GRINS-0.8.0
parsed_interior_qoi.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/multiphysics_sys.h"
31 #include "grins/assembly_context.h"
32 
33 // libMesh
34 #include "libmesh/getpot.h"
35 #include "libmesh/fem_system.h"
36 #include "libmesh/quadrature.h"
37 #include "libmesh/fe_base.h"
38 #include "libmesh/parsed_fem_function.h"
39 
40 namespace GRINS
41 {
42  ParsedInteriorQoI::ParsedInteriorQoI( const std::string& qoi_name )
43  : QoIBase(qoi_name) {}
44 
46  : QoIBase(original)
47  {
48  if (original.qoi_functional.get())
49  {
50  this->qoi_functional = original.qoi_functional->clone();
51  this->move_parameter
53  (original.qoi_functional.get()),
55  (this->qoi_functional.get()));
56  }
57  }
58 
60 
62  {
63  return new ParsedInteriorQoI( *this );
64  }
65 
67  (const GetPot& input,
68  const MultiphysicsSystem& system,
69  unsigned int /*qoi_num*/ )
70  {
73  (system, ""));
74  this->qoi_functional.reset(qf);
75 
76  this->set_parameter(*qf, input,
77  "QoI/ParsedInterior/qoi_functional", "DIE!");
78  }
79 
81  {
82  libMesh::FEBase* element_fe;
83  context.get_element_fe<libMesh::Real>(0, element_fe);
84  element_fe->get_JxW();
85  element_fe->get_xyz();
86 
87  qoi_functional->init_context(context);
88  }
89 
91  const unsigned int qoi_index )
92  {
93  libMesh::FEBase* element_fe;
94  context.get_element_fe<libMesh::Real>(0, element_fe);
95  const std::vector<libMesh::Real> &JxW = element_fe->get_JxW();
96 
97  const std::vector<libMesh::Point>& x_qp = element_fe->get_xyz();
98 
99  unsigned int n_qpoints = context.get_element_qrule().n_points();
100 
102  libMesh::Number& qoi = context.get_qois()[qoi_index];
103 
104  for( unsigned int qp = 0; qp != n_qpoints; qp++ )
105  {
106  const libMesh::Number func_val =
107  (*qoi_functional)(context, x_qp[qp], context.get_time());
108 
109  qoi += func_val * JxW[qp];
110  }
111  }
112 
114  const unsigned int qoi_index )
115  {
116  libMesh::FEBase* element_fe;
117  context.get_element_fe<libMesh::Real>(0, element_fe);
118  const std::vector<libMesh::Real> &JxW = element_fe->get_JxW();
119 
120  const std::vector<libMesh::Point>& x_qp = element_fe->get_xyz();
121 
122  // Local DOF count and quadrature point count
123  const unsigned int n_u_dofs = context.get_dof_indices().size();
124 
125  unsigned int n_qpoints = context.get_element_qrule().n_points();
126 
127  // Local solution vector - non-const version for finite
128  // differenting purposes
129  libMesh::DenseVector<libMesh::Number>& elem_solution =
130  const_cast<libMesh::DenseVector<libMesh::Number>&>
131  (context.get_elem_solution());
132 
134  libMesh::DenseVector<libMesh::Number> &Qu =
135  context.get_qoi_derivatives()[qoi_index];
136 
137  for( unsigned int qp = 0; qp != n_qpoints; qp++ )
138  {
139  // Central finite differencing to approximate derivatives.
140  // FIXME - we should hook the FParserAD stuff into
141  // ParsedFEMFunction
142 
143  for( unsigned int i = 0; i != n_u_dofs; ++i )
144  {
145  libMesh::Number &current_solution = elem_solution(i);
146  const libMesh::Number original_solution = current_solution;
147 
148  current_solution = original_solution + libMesh::TOLERANCE;
149 
150  const libMesh::Number plus_val =
151  (*qoi_functional)(context, x_qp[qp], context.get_time());
152 
153  current_solution = original_solution - libMesh::TOLERANCE;
154 
155  const libMesh::Number minus_val =
156  (*qoi_functional)(context, x_qp[qp], context.get_time());
157 
158  Qu(i) += (plus_val - minus_val) *
159  (0.5 / libMesh::TOLERANCE) * JxW[qp];
160 
161  // Don't forget to restore the correct solution...
162  current_solution = original_solution;
163  }
164  }
165  }
166 
167 } //namespace GRINS
ParsedInteriorQoI()
User never call default constructor.
GRINS namespace.
Parsed Interior QoI.
virtual void init(const GetPot &input, const MultiphysicsSystem &system, unsigned int qoi_num)
Initialize local variables.
Interface with libMesh for solving Multiphysics problems.
virtual void element_qoi(AssemblyContext &context, const unsigned int qoi_index)
Compute the qoi value.
virtual void element_qoi_derivative(AssemblyContext &context, const unsigned int qoi_index)
Compute the qoi derivative with respect to the solution.
virtual void init_context(AssemblyContext &context)
virtual QoIBase * clone() const
Required to provide clone (deep-copy) for adding QoI object to libMesh objects.
virtual void move_parameter(const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
When cloning an object, we need to update parameter pointers.
libMesh::UniquePtr< libMesh::FEMFunctionBase< libMesh::Number > > qoi_functional

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