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
40using namespace std;
41
42namespace Nektar
43{
46 "ImageWarpingSystem", ImageWarpingSystem::create,
47 "Image warping system.");
48
49/**
50 *
51 */
55 : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
56{
57}
58
59/**
60 *
61 */
62void ImageWarpingSystem::v_InitObject(bool DeclareField)
63{
65
66 // Define Velocity fields
68 int nq = m_fields[0]->GetNpoints();
69
70 for (int i = 0; i < m_spacedim; ++i)
71 {
73 }
74
75 // Bit of a hack: redefine u/v fields so they are continuous for
76 // Helmholtz solve.
79 m_session, m_graph, m_session->GetVariable(2));
80 m_fields[2] = fld;
82 *fld, m_graph, m_session->GetVariable(3));
83
84 // Tell UnsteadySystem to only integrate first two fields (i.e. I and
85 // phi).
86 m_intVariables.push_back(0);
87 m_intVariables.push_back(1);
88
89 // Define the normal velocity fields
90 if (m_fields[0]->GetTrace())
91 {
93 }
94
95 string advName;
96 string riemName;
97 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
100 m_advObject->SetFluxVector(&ImageWarpingSystem::GetFluxVector, this);
101 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
103 riemName, m_session);
105 this);
106
107 m_advObject->SetRiemannSolver(m_riemannSolver);
108 m_advObject->InitObject(m_session, m_fields);
109
111 {
114 }
115 else
116 {
117 ASSERTL0(false, "Implicit unsteady Advection not set up.");
118 }
119}
120
121/**
122 *
123 */
125{
126}
127
128/**
129 *
130 */
132 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
133 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
134{
135 boost::ignore_unused(time);
136
137 int npoints = GetNpoints();
138 int ncoeffs = inarray[0].size();
140
141 // Load parameter alpha.
142 m_session->LoadParameter("Alpha", m_alpha);
143
145 "CG not implemented yet.");
146
147 // Set up storage arrays.
148 Array<OneD, NekDouble> tmp(npoints);
149 Array<OneD, NekDouble> alloc(3 * npoints);
150 Array<OneD, NekDouble> dIdx1(alloc);
151 Array<OneD, NekDouble> dIdx2(alloc + npoints);
152 Array<OneD, NekDouble> dIdx3(alloc + 2 * npoints);
153
154 // Calculate grad I.
155 m_fields[0]->PhysDeriv(inarray[0], dIdx1, dIdx2);
156
157 // Set factors.
158 // TODO: Check - should be -1?
160
161 // Multiply by phi, and perform Helmholtz solve to calculate the
162 // advection velocity field.
163 for (int i = 0; i < 2; ++i)
164 {
165 Vmath::Vmul(npoints, &alloc[i * npoints], 1, inarray[1].get(), 1,
166 m_fields[i + 2]->UpdatePhys().get(), 1);
167 Vmath::Smul(npoints, 1 / m_alpha / m_alpha,
168 m_fields[i + 2]->GetPhys().get(), 1,
169 m_fields[i + 2]->UpdatePhys().get(), 1);
170 m_fields[i + 2]->HelmSolve(m_fields[i + 2]->GetPhys(),
171 m_fields[i + 2]->UpdateCoeffs(), factors);
172 m_fields[i + 2]->BwdTrans(m_fields[i + 2]->GetCoeffs(), m_velocity[i]);
173 }
174
175 // Calculate the weak advection operator for I and phi - result is put
176 // in WeakAdv and is in physical space.
177 m_advObject->Advect(2, m_fields, m_velocity, inarray, outarray, 0.0);
178 for (int i = 0; i < 2; ++i)
179 {
180 Vmath::Neg(npoints, outarray[i], 1);
181 }
182
183 // Calculate du/dx -> dIdx1, dv/dy -> dIdx2.
184 m_fields[2]->PhysDeriv(m_velocity[0], dIdx1, dIdx3);
185 m_fields[3]->PhysDeriv(m_velocity[1], dIdx3, dIdx2);
186
187 // Calculate RHS = I*div(u) = I*du/dx + I*dv/dy -> dIdx1.
188 Vmath::Vvtvvtp(npoints, dIdx1.get(), 1, inarray[0].get(), 1, dIdx2.get(), 1,
189 inarray[0].get(), 1, dIdx1.get(), 1);
190
191 // Take inner product to get to coefficient space.
192 Array<OneD, NekDouble> tmp2(ncoeffs);
193 m_fields[0]->IProductWRTBase(dIdx1, tmp2);
194
195 // Multiply by elemental inverse mass matrix, backwards transform
196 m_fields[0]->MultiplyByElmtInvMass(tmp2, tmp2);
197 m_fields[0]->BwdTrans(tmp2, tmp);
198 Vmath::Vadd(npoints, outarray[0], 1, tmp, 1, outarray[0], 1);
199}
200
201/**
202 *
203 */
205 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
206 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
207{
208 int nvariables = inarray.size();
210
211 switch (m_projectionType)
212 {
214 {
215 // Just copy over array
216 if (inarray != outarray)
217 {
218 int npoints = GetNpoints();
219
220 for (int i = 0; i < nvariables; ++i)
221 {
222 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
223 }
224 }
225 }
226 break;
227 default:
228 ASSERTL0(false, "Unknown projection scheme");
229 break;
230 }
231}
232
233/**
234 * @brief Get the normal velocity
235 */
237{
238 // Number of trace (interface) points
239 int nTracePts = GetTraceNpoints();
240
241 // Auxiliary variable to compute the normal velocity
242 Array<OneD, NekDouble> tmp(nTracePts);
243
244 // Reset the normal velocity
245 Vmath::Zero(nTracePts, m_traceVn, 1);
246
247 for (int i = 0; i < m_velocity.size(); ++i)
248 {
249 m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
250
251 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
252 m_traceVn, 1);
253 }
254
255 return m_traceVn;
256}
257
258/**
259 *
260 */
262 const Array<OneD, Array<OneD, NekDouble>> &physfield,
264{
265 ASSERTL1(flux[0].size() == m_velocity.size(),
266 "Dimension of flux array and velocity array do not match");
267
268 int nq = physfield[0].size();
269
270 for (int i = 0; i < flux.size(); ++i)
271 {
272 for (int j = 0; j < flux[0].size(); ++j)
273 {
274 Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
275 }
276 }
277}
278
279/**
280 *
281 */
283{
285}
286} // 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
ImageWarpingSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Array< OneD, Array< OneD, NekDouble > > m_velocity
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
static std::string className
Name of class.
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
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
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:270
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:176
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:408
StdRegions::ConstFactorMap factors
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:207
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513
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:569
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:354
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:245
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
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:687
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191