MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
system_initialization.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2020 Manav Bhatia and MAST authors
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 // MAST includes
22 #include "base/nonlinear_system.h"
23 #include "numerics/utility.h"
25 
26 // libMesh includes
27 #include "libmesh/function_base.h"
28 
29 
30 
31 
33  const std::string& prefix):
34 _system(sys),
35 _prefix(prefix) {
36 
37  // initialize the point locator for this mesh
38  sys.system().get_mesh().sub_point_locator();
39 }
40 
41 
42 
44 { }
45 
46 
47 
48 
49 unsigned int
51  return _system.n_vars();
52 }
53 
54 
55 
56 
57 const libMesh::FEType&
58 MAST::SystemInitialization::fetype(unsigned int i) const {
59 
60  return _system.variable_type(i);
61 }
62 
63 
64 
65 void
67 
68  // make sure that the dimension of the sol vector matches the dimension
69  // specified for this system
70  libmesh_assert_equal_to(sol.size(), _vars.size());
71 
72  // now create a function and use it for initialization
73  class SolutionFunction:
74  public libMesh::FunctionBase<Real> {
75  public:
76  SolutionFunction(const RealVectorX& s):
77  libMesh::FunctionBase<Real>() {
78  _sol.resize((unsigned int)s.size());
79  for (unsigned int i=0; i<s.size(); i++) _sol(i) = s(i);
80  }
81 
82 
83  SolutionFunction(const DenseRealVector& sol):
84  libMesh::FunctionBase<Real>(),
85  _sol(sol) { }
86 
87  virtual std::unique_ptr<libMesh::FunctionBase<Real> > clone () const {
88  libMesh::FunctionBase<Real> *rval = new SolutionFunction(_sol);
89  return std::unique_ptr<libMesh::FunctionBase<Real> >(rval);
90  }
91 
92  // this should not get called
93  virtual Real operator()
94  (const libMesh::Point& p, const Real time) {libmesh_assert(false);}
95 
96  virtual void
97  operator() (const libMesh::Point& p,
98  const Real time,
99  libMesh::DenseVector<Real>& output) {
100  output = _sol;
101  }
102  protected:
103  DenseRealVector _sol;
104  };
105 
106  SolutionFunction sol_func(sol);
107 
108  _system.project_solution(&sol_func);
109 }
110 
111 
112 
113 void
116 
117  // now create a function and use it for initialization
118  class SolutionFunction:
119  public libMesh::FunctionBase<Real> {
120  public:
121  SolutionFunction(const unsigned int n_vars,
123  libMesh::FunctionBase<Real>(),
124  _nvars(n_vars),
125  _sol(s) { }
126 
127  virtual std::unique_ptr<libMesh::FunctionBase<Real> > clone () const {
128  libMesh::FunctionBase<Real> *rval = new SolutionFunction(_nvars, _sol);
129  return std::unique_ptr<libMesh::FunctionBase<Real> >(rval);
130  }
131 
132  // this should not get called
133  virtual Real operator()
134  (const libMesh::Point& p, const Real time) {libmesh_assert(false);}
135 
136  virtual void
137  operator() (const libMesh::Point& p,
138  const Real time,
139  libMesh::DenseVector<Real>& output) {
140 
141  RealVectorX v = RealVectorX::Zero(_nvars);
142  _sol(p, time, v);
143  MAST::copy(output, v);
144  }
145  protected:
146  const unsigned int _nvars;
148  };
149 
150  SolutionFunction sol_func(_vars.size(), sol);
151 
152  _system.project_solution(&sol_func);
153 }
154 
155 
void initialize_solution(const RealVectorX &sol)
initializes the FE solution vector to the constant solution provided in sol.
std::vector< unsigned int > _vars
This class implements a system for solution of nonlinear systems.
virtual ~SystemInitialization()
virtual destructor
const libMesh::FEType & fetype(unsigned int i) const
libMesh::Real Real
SystemInitialization(MAST::NonlinearSystem &sys, const std::string &prefix)
initialize the variables in the provided system sys of order and family.
libMesh::DenseVector< Real > DenseRealVector
void copy(DenseRealMatrix &m1, const RealMatrixX &m2)
Definition: utility.h:167
Matrix< Real, Dynamic, 1 > RealVectorX
MAST::NonlinearSystem & _system