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 
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44  string ImageWarpingSystem::className =
46  "ImageWarpingSystem",
47  ImageWarpingSystem::create,
48  "Image warping system.");
49 
50  ImageWarpingSystem::ImageWarpingSystem(
53  : UnsteadySystem(pSession, pGraph),
54  AdvectionSystem(pSession, pGraph)
55  {
56  }
57 
59  {
61 
62  // Define Velocity fields
64  int nq = m_fields[0]->GetNpoints();
65 
66  for(int i = 0; i < m_spacedim; ++i)
67  {
68  m_velocity[i] = Array<OneD, NekDouble> (nq,0.0);
69  }
70 
71  // Bit of a hack: redefine u/v fields so they are continuous for
72  // Helmholtz solve.
76  m_fields[2] = fld;
78  ::AllocateSharedPtr(*fld,m_graph,m_session->GetVariable(3));
79 
80  // Tell UnsteadySystem to only integrate first two fields (i.e. I and
81  // phi).
82  m_intVariables.push_back(0);
83  m_intVariables.push_back(1);
84 
85  // Define the normal velocity fields
86  if (m_fields[0]->GetTrace())
87  {
89  }
90 
91  string advName;
92  string riemName;
93  m_session->LoadSolverInfo(
94  "AdvectionType", advName, "WeakDG");
96  GetAdvectionFactory().CreateInstance(advName, advName);
97  m_advObject->SetFluxVector(
99  m_session->LoadSolverInfo(
100  "UpwindType", riemName, "Upwind");
103  m_riemannSolver->SetScalar(
105 
106  m_advObject->SetRiemannSolver(m_riemannSolver);
107  m_advObject->InitObject(m_session, m_fields);
108 
110  {
113  }
114  else
115  {
116  ASSERTL0(false, "Implicit unsteady Advection not set up.");
117  }
118  }
119 
121  {
122 
123  }
124 
126  const Array<OneD, const Array<OneD,NekDouble> > &inarray,
127  Array<OneD, Array<OneD,NekDouble> > &outarray,
128  const NekDouble time)
129  {
130  boost::ignore_unused(time);
131 
132  int npoints = GetNpoints();
133  int ncoeffs = inarray[0].num_elements();
135 
136  // Load parameter alpha.
137  m_session->LoadParameter("Alpha", m_alpha);
138 
140  "CG not implemented yet.");
141 
142  // Set up storage arrays.
143  Array<OneD, NekDouble> tmp(npoints);
144  Array<OneD, NekDouble> alloc(3*npoints);
145  Array<OneD, NekDouble> dIdx1(alloc);
146  Array<OneD, NekDouble> dIdx2(alloc+npoints);
147  Array<OneD, NekDouble> dIdx3(alloc+2*npoints);
148 
149  // Calculate grad I.
150  m_fields[0]->PhysDeriv(inarray[0], dIdx1, dIdx2);
151 
152  // Set factors.
153  // TODO: Check - should be -1?
154  factors[StdRegions::eFactorLambda] = 1.0 / m_alpha / m_alpha;
155 
156  // Multiply by phi, and perform Helmholtz solve to calculate the
157  // advection velocity field.
158  for (int i = 0; i < 2; ++i)
159  {
160  Vmath::Vmul(npoints, &alloc[i*npoints], 1, inarray[1].get(), 1,
161  m_fields[i+2]->UpdatePhys().get(), 1);
162  Vmath::Smul(npoints, 1/m_alpha/m_alpha, m_fields[i+2]->GetPhys().get(), 1,
163  m_fields[i+2]->UpdatePhys().get(), 1);
164  m_fields[i+2]->HelmSolve(m_fields[i+2]->GetPhys(),
165  m_fields[i+2]->UpdateCoeffs(),
166  NullFlagList, factors);
167  m_fields[i+2]->BwdTrans(m_fields[i+2]->GetCoeffs(),
168  m_velocity[i]);
169  }
170 
171  // Calculate the weak advection operator for I and phi - result is put
172  // in WeakAdv and is in physical space.
173  m_advObject->Advect(2, m_fields, m_velocity, inarray,
174  outarray, 0.0);
175  for(int i = 0; i < 2; ++i)
176  {
177  Vmath::Neg(npoints, outarray[i], 1);
178  }
179 
180  // Calculate du/dx -> dIdx1, dv/dy -> dIdx2.
181  m_fields[2]->PhysDeriv(m_velocity[0], dIdx1, dIdx3);
182  m_fields[3]->PhysDeriv(m_velocity[1], dIdx3, dIdx2);
183 
184  // Calculate RHS = I*div(u) = I*du/dx + I*dv/dy -> dIdx1.
185  Vmath::Vvtvvtp(npoints, dIdx1.get(), 1, inarray[0].get(), 1,
186  dIdx2.get(), 1, inarray[0].get(), 1,
187  dIdx1.get(), 1);
188 
189  // Take inner product to get to coefficient space.
190  Array<OneD, NekDouble> tmp2(ncoeffs);
191  m_fields[0]->IProductWRTBase(dIdx1, tmp2);
192 
193  // Multiply by elemental inverse mass matrix, backwards transform
194  m_fields[0]->MultiplyByElmtInvMass(tmp2, tmp2);
195  m_fields[0]->BwdTrans(tmp2,tmp);
196  Vmath::Vadd(npoints, outarray[0], 1, tmp, 1, outarray[0], 1);
197  }
198 
199 
200 
201  /**
202  *
203  */
205  const Array<OneD, NekDouble> >&inarray,
206  Array<OneD, Array<OneD, NekDouble> >&outarray,
207  const NekDouble time)
208  {
209  int nvariables = inarray.num_elements();
210  SetBoundaryConditions(time);
211 
212  switch(m_projectionType)
213  {
215  {
216  // Just copy over array
217  int npoints = GetNpoints();
218 
219  for(int i = 0; i < nvariables; ++i)
220  {
221  Vmath::Vcopy(npoints,inarray[i],1,outarray[i],1);
222  }
223  }
224  break;
225  default:
226  ASSERTL0(false,"Unknown projection scheme");
227  break;
228  }
229  }
230 
231  /**
232  * @brief Get the normal velocity
233  */
235  {
236  // Number of trace (interface) points
237  int nTracePts = GetTraceNpoints();
238 
239  // Auxiliary variable to compute the normal velocity
240  Array<OneD, NekDouble> tmp(nTracePts);
241 
242  // Reset the normal velocity
243  Vmath::Zero(nTracePts, m_traceVn, 1);
244 
245  for (int i = 0; i < m_velocity.num_elements(); ++i)
246  {
247  m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
248 
249  Vmath::Vvtvp(nTracePts,
250  m_traceNormals[i], 1,
251  tmp, 1,
252  m_traceVn, 1,
253  m_traceVn, 1);
254  }
255 
256  return m_traceVn;
257  }
258 
260  const Array<OneD, Array<OneD, NekDouble> > &physfield,
262  {
263  ASSERTL1(flux[0].num_elements() == m_velocity.num_elements(),
264  "Dimension of flux array and velocity array do not match");
265 
266  int nq = physfield[0].num_elements();
267 
268  for (int i = 0; i < flux.num_elements(); ++i)
269  {
270  for (int j = 0; j < flux[0].num_elements(); ++j)
271  {
272  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1,
273  flux[i][j], 1);
274  }
275  }
276  }
277 
279  {
281  }
282 }
Array< OneD, Array< OneD, NekDouble > > m_velocity
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:163
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
std::shared_ptr< ContField2D > ContField2DSharedPtr
Definition: ContField2D.h:289
virtual void v_InitObject()
Init object for UnsteadySystem class.
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:445
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
STL namespace.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
Array< OneD, NekDouble > m_traceVn
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Print a summary of time stepping parameters.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
Base class for unsteady solvers.
RiemannSolverFactory & GetRiemannSolverFactory()
int m_spacedim
Spatial dimension (>= expansion dim).
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
double NekDouble
EquationSystemFactory & GetEquationSystemFactory()
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT int GetNpoints()
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:540
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
A base class for PDEs which include an advection component.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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:302
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:186
static FlagList NullFlagList
An empty flag list.