www.mooseframework.org
MooseVariableInterface.C
Go to the documentation of this file.
1 /****************************************************************/
2 /* DO NOT MODIFY THIS HEADER */
3 /* MOOSE - Multiphysics Object Oriented Simulation Environment */
4 /* */
5 /* (c) 2010 Battelle Energy Alliance, LLC */
6 /* ALL RIGHTS RESERVED */
7 /* */
8 /* Prepared by Battelle Energy Alliance, LLC */
9 /* Under Contract No. DE-AC07-05ID14517 */
10 /* With the U. S. Department of Energy */
11 /* */
12 /* See COPYRIGHT for full restrictions */
13 /****************************************************************/
14 
15 #include "MooseVariableInterface.h"
16 
17 #include "Assembly.h"
18 #include "MooseError.h" // mooseDeprecated
19 #include "MooseTypes.h"
20 #include "MooseVariable.h"
21 #include "Problem.h"
22 #include "SubProblem.h"
23 
25  bool nodal,
26  std::string var_param_name)
27  : _nodal(nodal)
28 {
29  const InputParameters & parameters = moose_object->parameters();
30 
31  SubProblem & problem = *parameters.get<SubProblem *>("_subproblem");
32 
33  THREAD_ID tid = parameters.get<THREAD_ID>("_tid");
34 
35  // Try the scalar version first
36  std::string variable_name = parameters.getMooseType(var_param_name);
37  if (variable_name == "")
38  // When using vector variables, we are only going to use the first one in the list at the
39  // interface level...
40  variable_name = parameters.getVecMooseType(var_param_name)[0];
41 
42  _variable = &problem.getVariable(tid, variable_name);
43 
44  _mvi_assembly = &problem.assembly(tid);
45 }
46 
48 
51 {
52  return _variable;
53 }
54 
55 const VariableValue &
57 {
58  if (_nodal)
59  return _variable->nodalSln();
60  else
61  return _variable->sln();
62 }
63 
64 const VariableValue &
66 {
67  if (_nodal)
68  return _variable->nodalSlnOld();
69  else
70  return _variable->slnOld();
71 }
72 
73 const VariableValue &
75 {
76  if (_nodal)
77  return _variable->nodalSlnOlder();
78  else
79  return _variable->slnOlder();
80 }
81 
82 const VariableValue &
84 {
85  if (_nodal)
86  return _variable->nodalSlnDot();
87  else
88  return _variable->uDot();
89 }
90 
91 const VariableValue &
93 {
94  if (_nodal)
95  return _variable->nodalSlnDuDotDu();
96  else
97  return _variable->duDotDu();
98 }
99 
100 const VariableGradient &
102 {
103  if (_nodal)
104  mooseError("Nodal variables do not have gradients");
105 
106  return _variable->gradSln();
107 }
108 
109 const VariableGradient &
111 {
112  if (_nodal)
113  mooseError("Nodal variables do not have gradients");
114 
115  return _variable->gradSlnOld();
116 }
117 
118 const VariableGradient &
120 {
121  if (_nodal)
122  mooseError("Nodal variables do not have gradients");
123 
124  return _variable->gradSlnOlder();
125 }
126 
127 const VariableSecond &
129 {
130  if (_nodal)
131  mooseError("Nodal variables do not have second derivatives");
132 
133  return _variable->secondSln();
134 }
135 
136 const VariableSecond &
138 {
139  if (_nodal)
140  mooseError("Nodal variables do not have second derivatives");
141 
142  return _variable->secondSlnOld();
143 }
144 
145 const VariableSecond &
147 {
148  if (_nodal)
149  mooseError("Nodal variables do not have second derivatives");
150 
151  return _variable->secondSlnOlder();
152 }
153 
154 const VariableTestSecond &
156 {
157  if (_nodal)
158  mooseError("Nodal variables do not have second derivatives");
159 
160  return _variable->secondPhi();
161 }
162 
163 const VariableTestSecond &
165 {
166  if (_nodal)
167  mooseError("Nodal variables do not have second derivatives");
168 
169  return _variable->secondPhiFace();
170 }
171 
172 const VariablePhiSecond &
174 {
175  if (_nodal)
176  mooseError("Nodal variables do not have second derivatives");
177 
178  return _mvi_assembly->secondPhi();
179 }
180 
181 const VariablePhiSecond &
183 {
184  if (_nodal)
185  mooseError("Nodal variables do not have second derivatives");
186 
187  return _mvi_assembly->secondPhiFace();
188 }
const VariablePhiSecond & secondPhiFace()
const VariableGradient & gradSlnOld()
virtual const VariableTestSecond & secondTestFace()
The second derivative of the test function on the current face.
virtual const VariableSecond & secondOlder()
The older second derivative of the variable this object is operating on.
virtual const VariableValue & valueOld()
The old value of the variable this object is operating on.
Class for stuff related to variables.
Definition: MooseVariable.h:43
const VariableValue & slnOlder()
virtual const VariableGradient & gradientOld()
The old gradient of the variable this object is operating on.
virtual MooseVariable & getVariable(THREAD_ID tid, const std::string &var_name)=0
Returns the variable reference for requested variable which may be in any system. ...
std::vector< std::string > getVecMooseType(const std::string &name) const
void mooseError(Args &&...args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:182
const VariableSecond & secondSln()
virtual Assembly & assembly(THREAD_ID tid)=0
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
bool _nodal
Whether or not this object is acting only at nodes.
const VariableSecond & secondSlnOld()
const VariableSecond & secondSlnOlder()
virtual const VariableGradient & gradientOlder()
The older gradient of the variable this object is operating on.
const VariablePhiSecond & secondPhi()
Definition: Assembly.h:536
virtual const VariableSecond & second()
The second derivative of the variable this object is operating on.
const VariableValue & nodalSlnDot()
virtual const VariableValue & dot()
The time derivative of the variable this object is operating on.
MooseVariableInterface(const MooseObject *moose_object, bool nodal, std::string var_param_name="variable")
Constructing the object.
MooseVariable * mooseVariable()
Get the variable that this object is using.
virtual const VariableValue & valueOlder()
The older value of the variable this object is operating on.
virtual const VariablePhiSecond & secondPhi()
The second derivative of the trial function.
const VariablePhiSecond & secondPhi()
const VariableValue & slnOld()
Every object that can be built by the factory should be derived from this class.
Definition: MooseObject.h:36
const VariableValue & nodalSlnOld()
const VariableValue & duDotDu()
virtual const VariablePhiSecond & secondPhiFace()
The second derivative of the trial function on the current face.
std::string getMooseType(const std::string &name) const
Utility functions for retrieving one of the MooseTypes variables into the common "string" base class...
const VariableValue & nodalSlnOlder()
const VariableValue & nodalSln()
const VariableValue & nodalSlnDuDotDu()
Generic class for solving transient nonlinear problems.
Definition: SubProblem.h:53
virtual const VariableValue & value()
The value of the variable this object is operating on.
const VariableValue & sln()
const VariableValue & uDot()
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:53
const VariableGradient & gradSlnOlder()
const VariablePhiSecond & secondPhiFace()
Definition: Assembly.h:540
virtual const VariableValue & dotDu()
The derivative of the time derivative of the variable this object is operating on with respect to thi...
virtual const VariableSecond & secondOld()
The old second derivative of the variable this object is operating on.
MooseVariable * _variable
The variable this object is acting on.
virtual const VariableGradient & gradient()
The gradient of the variable this object is operating on.
const VariableGradient & gradSln()
virtual const VariableTestSecond & secondTest()
The second derivative of the test function.
unsigned int THREAD_ID
Definition: MooseTypes.h:79