Nektar++
ImageWarpingSystem.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ImageWarpingSystem.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Image warping solve routines
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
38 #include <MultiRegions/ContField.h>
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44 string ImageWarpingSystem::className =
46  "ImageWarpingSystem", ImageWarpingSystem::create,
47  "Image warping system.");
48 
49 ImageWarpingSystem::ImageWarpingSystem(
52  : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
53 {
54 }
55 
56 void ImageWarpingSystem::v_InitObject(bool DeclareField)
57 {
58  AdvectionSystem::v_InitObject(DeclareField);
59 
60  // Define Velocity fields
62  int nq = m_fields[0]->GetNpoints();
63 
64  for (int i = 0; i < m_spacedim; ++i)
65  {
66  m_velocity[i] = Array<OneD, NekDouble>(nq, 0.0);
67  }
68 
69  // Bit of a hack: redefine u/v fields so they are continuous for
70  // Helmholtz solve.
73  m_session, m_graph, m_session->GetVariable(2));
74  m_fields[2] = fld;
76  *fld, m_graph, m_session->GetVariable(3));
77 
78  // Tell UnsteadySystem to only integrate first two fields (i.e. I and
79  // phi).
80  m_intVariables.push_back(0);
81  m_intVariables.push_back(1);
82 
83  // Define the normal velocity fields
84  if (m_fields[0]->GetTrace())
85  {
87  }
88 
89  string advName;
90  string riemName;
91  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
92  m_advObject =
94  m_advObject->SetFluxVector(&ImageWarpingSystem::GetFluxVector, this);
95  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
97  riemName, m_session);
99  this);
100 
101  m_advObject->SetRiemannSolver(m_riemannSolver);
102  m_advObject->InitObject(m_session, m_fields);
103 
105  {
108  }
109  else
110  {
111  ASSERTL0(false, "Implicit unsteady Advection not set up.");
112  }
113 }
114 
116 {
117 }
118 
120  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
121  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
122 {
123  boost::ignore_unused(time);
124 
125  int npoints = GetNpoints();
126  int ncoeffs = inarray[0].size();
128 
129  // Load parameter alpha.
130  m_session->LoadParameter("Alpha", m_alpha);
131 
133  "CG not implemented yet.");
134 
135  // Set up storage arrays.
136  Array<OneD, NekDouble> tmp(npoints);
137  Array<OneD, NekDouble> alloc(3 * npoints);
138  Array<OneD, NekDouble> dIdx1(alloc);
139  Array<OneD, NekDouble> dIdx2(alloc + npoints);
140  Array<OneD, NekDouble> dIdx3(alloc + 2 * npoints);
141 
142  // Calculate grad I.
143  m_fields[0]->PhysDeriv(inarray[0], dIdx1, dIdx2);
144 
145  // Set factors.
146  // TODO: Check - should be -1?
147  factors[StdRegions::eFactorLambda] = 1.0 / m_alpha / m_alpha;
148 
149  // Multiply by phi, and perform Helmholtz solve to calculate the
150  // advection velocity field.
151  for (int i = 0; i < 2; ++i)
152  {
153  Vmath::Vmul(npoints, &alloc[i * npoints], 1, inarray[1].get(), 1,
154  m_fields[i + 2]->UpdatePhys().get(), 1);
155  Vmath::Smul(npoints, 1 / m_alpha / m_alpha,
156  m_fields[i + 2]->GetPhys().get(), 1,
157  m_fields[i + 2]->UpdatePhys().get(), 1);
158  m_fields[i + 2]->HelmSolve(m_fields[i + 2]->GetPhys(),
159  m_fields[i + 2]->UpdateCoeffs(), factors);
160  m_fields[i + 2]->BwdTrans(m_fields[i + 2]->GetCoeffs(), m_velocity[i]);
161  }
162 
163  // Calculate the weak advection operator for I and phi - result is put
164  // in WeakAdv and is in physical space.
165  m_advObject->Advect(2, m_fields, m_velocity, inarray, outarray, 0.0);
166  for (int i = 0; i < 2; ++i)
167  {
168  Vmath::Neg(npoints, outarray[i], 1);
169  }
170 
171  // Calculate du/dx -> dIdx1, dv/dy -> dIdx2.
172  m_fields[2]->PhysDeriv(m_velocity[0], dIdx1, dIdx3);
173  m_fields[3]->PhysDeriv(m_velocity[1], dIdx3, dIdx2);
174 
175  // Calculate RHS = I*div(u) = I*du/dx + I*dv/dy -> dIdx1.
176  Vmath::Vvtvvtp(npoints, dIdx1.get(), 1, inarray[0].get(), 1, dIdx2.get(), 1,
177  inarray[0].get(), 1, dIdx1.get(), 1);
178 
179  // Take inner product to get to coefficient space.
180  Array<OneD, NekDouble> tmp2(ncoeffs);
181  m_fields[0]->IProductWRTBase(dIdx1, tmp2);
182 
183  // Multiply by elemental inverse mass matrix, backwards transform
184  m_fields[0]->MultiplyByElmtInvMass(tmp2, tmp2);
185  m_fields[0]->BwdTrans(tmp2, tmp);
186  Vmath::Vadd(npoints, outarray[0], 1, tmp, 1, outarray[0], 1);
187 }
188 
189 /**
190  *
191  */
193  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
194  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
195 {
196  int nvariables = inarray.size();
197  SetBoundaryConditions(time);
198 
199  switch (m_projectionType)
200  {
202  {
203  // Just copy over array
204  int npoints = GetNpoints();
205 
206  for (int i = 0; i < nvariables; ++i)
207  {
208  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
209  }
210  }
211  break;
212  default:
213  ASSERTL0(false, "Unknown projection scheme");
214  break;
215  }
216 }
217 
218 /**
219  * @brief Get the normal velocity
220  */
222 {
223  // Number of trace (interface) points
224  int nTracePts = GetTraceNpoints();
225 
226  // Auxiliary variable to compute the normal velocity
227  Array<OneD, NekDouble> tmp(nTracePts);
228 
229  // Reset the normal velocity
230  Vmath::Zero(nTracePts, m_traceVn, 1);
231 
232  for (int i = 0; i < m_velocity.size(); ++i)
233  {
234  m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
235 
236  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
237  m_traceVn, 1);
238  }
239 
240  return m_traceVn;
241 }
242 
244  const Array<OneD, Array<OneD, NekDouble>> &physfield,
246 {
247  ASSERTL1(flux[0].size() == m_velocity.size(),
248  "Dimension of flux array and velocity array do not match");
249 
250  int nq = physfield[0].size();
251 
252  for (int i = 0; i < flux.size(); ++i)
253  {
254  for (int j = 0; j < flux[0].size(); ++j)
255  {
256  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
257  }
258  }
259 }
260 
262 {
264 }
265 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
Array< OneD, Array< OneD, NekDouble > > m_velocity
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &flux)
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
virtual void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
Array< OneD, NekDouble > m_traceVn
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
A base class for PDEs which include an advection component.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT int GetNpoints()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:289
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:48
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:574
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:359
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:248
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.cpp:692
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255