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
37
38using namespace std;
39
40namespace Nektar
41{
44 "ImageWarpingSystem", ImageWarpingSystem::create,
45 "Image warping system.");
46
47/**
48 *
49 */
53 : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
54{
55}
56
57/**
58 *
59 */
60void ImageWarpingSystem::v_InitObject(bool DeclareField)
61{
63
64 // Define Velocity fields
66 int nq = m_fields[0]->GetNpoints();
67
68 for (int i = 0; i < m_spacedim; ++i)
69 {
71 }
72
73 // Bit of a hack: redefine u/v fields so they are continuous for
74 // Helmholtz solve.
77 m_session, m_graph, m_session->GetVariable(2));
78 m_fields[2] = fld;
80 *fld, m_graph, m_session->GetVariable(3));
81
82 // Tell UnsteadySystem to only integrate first two fields (i.e. I and
83 // phi).
84 m_intVariables.push_back(0);
85 m_intVariables.push_back(1);
86
87 // Define the normal velocity fields
88 if (m_fields[0]->GetTrace())
89 {
91 }
92
93 string advName;
94 string riemName;
95 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
99 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
101 riemName, m_session);
103 this);
104
105 m_advObject->SetRiemannSolver(m_riemannSolver);
106 m_advObject->InitObject(m_session, m_fields);
107
109 {
112 }
113 else
114 {
115 ASSERTL0(false, "Implicit unsteady Advection not set up.");
116 }
117}
118
119/**
120 *
121 */
123{
124}
125
126/**
127 *
128 */
130 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
132 [[maybe_unused]] const NekDouble time)
133{
134 int npoints = GetNpoints();
135 int ncoeffs = inarray[0].size();
137
138 // Load parameter alpha.
139 m_session->LoadParameter("Alpha", m_alpha);
140
142 "CG not implemented yet.");
143
144 // Set up storage arrays.
145 Array<OneD, NekDouble> tmp(npoints);
146 Array<OneD, NekDouble> alloc(3 * npoints);
147 Array<OneD, NekDouble> dIdx1(alloc);
148 Array<OneD, NekDouble> dIdx2(alloc + npoints);
149 Array<OneD, NekDouble> dIdx3(alloc + 2 * npoints);
150
151 // Calculate grad I.
152 m_fields[0]->PhysDeriv(inarray[0], dIdx1, dIdx2);
153
154 // Set factors.
155 // TODO: Check - should be -1?
157
158 // Multiply by phi, and perform Helmholtz solve to calculate the
159 // advection velocity field.
160 for (int i = 0; i < 2; ++i)
161 {
162 Vmath::Vmul(npoints, &alloc[i * npoints], 1, inarray[1].get(), 1,
163 m_fields[i + 2]->UpdatePhys().get(), 1);
164 Vmath::Smul(npoints, 1 / m_alpha / m_alpha,
165 m_fields[i + 2]->GetPhys().get(), 1,
166 m_fields[i + 2]->UpdatePhys().get(), 1);
167 m_fields[i + 2]->HelmSolve(m_fields[i + 2]->GetPhys(),
168 m_fields[i + 2]->UpdateCoeffs(), factors);
169 m_fields[i + 2]->BwdTrans(m_fields[i + 2]->GetCoeffs(), m_velocity[i]);
170 }
171
172 // Calculate the weak advection operator for I and phi - result is put
173 // in WeakAdv and is in physical space.
174 m_advObject->Advect(2, m_fields, m_velocity, inarray, 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, dIdx2.get(), 1,
186 inarray[0].get(), 1, dIdx1.get(), 1);
187
188 // Take inner product to get to coefficient space.
189 Array<OneD, NekDouble> tmp2(ncoeffs);
190 m_fields[0]->IProductWRTBase(dIdx1, tmp2);
191
192 // Multiply by elemental inverse mass matrix, backwards transform
193 m_fields[0]->MultiplyByElmtInvMass(tmp2, tmp2);
194 m_fields[0]->BwdTrans(tmp2, tmp);
195 Vmath::Vadd(npoints, outarray[0], 1, tmp, 1, outarray[0], 1);
196}
197
198/**
199 *
200 */
202 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
203 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
204{
205 int nvariables = inarray.size();
207
208 switch (m_projectionType)
209 {
211 {
212 // Just copy over array
213 if (inarray != outarray)
214 {
215 int npoints = GetNpoints();
216
217 for (int i = 0; i < nvariables; ++i)
218 {
219 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
220 }
221 }
222 }
223 break;
224 default:
225 ASSERTL0(false, "Unknown projection scheme");
226 break;
227 }
228}
229
230/**
231 * @brief Get the normal velocity
232 */
234{
235 // Number of trace (interface) points
236 int nTracePts = GetTraceNpoints();
237
238 // Auxiliary variable to compute the normal velocity
239 Array<OneD, NekDouble> tmp(nTracePts);
240
241 // Reset the normal velocity
242 Vmath::Zero(nTracePts, m_traceVn, 1);
243
244 for (int i = 0; i < m_velocity.size(); ++i)
245 {
246 m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
247
248 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
249 m_traceVn, 1);
250 }
251
252 return m_traceVn;
253}
254
255/**
256 *
257 */
259 const Array<OneD, Array<OneD, NekDouble>> &physfield,
261{
262 ASSERTL1(flux[0].size() == m_velocity.size(),
263 "Dimension of flux array and velocity array do not match");
264
265 int nq = physfield[0].size();
266
267 for (int i = 0; i < flux.size(); ++i)
268 {
269 for (int j = 0; j < flux[0].size(); ++j)
270 {
271 Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
272 }
273 }
274}
275
276/**
277 *
278 */
280{
282}
283} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
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)
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.
void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
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.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
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.
SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
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.
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:268
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
StdRegions::ConstFactorMap factors
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.hpp:72
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
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.hpp:366
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.hpp:180
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.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
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.hpp:439
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.