Nektar++
SubSteppingExtrapolate.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: SubSteppingExtrapolate.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: Abstract base class for SubSteppingExtrapolate.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
36 
38 
39 namespace Nektar
40 {
41 using namespace LibUtilities;
42 
43 /**
44  * Registers the class with the Factory.
45  */
48  "SubStepping", SubSteppingExtrapolate::create, "SubStepping");
49 
54  const SolverUtils::AdvectionSharedPtr advObject)
55  : Extrapolate(pSession, pFields, pPressure, pVel, advObject)
56 {
57  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
58  m_session->LoadParameter("SubStepCFL", m_cflSafetyFactor, 0.5);
59  m_session->LoadParameter("MinSubSteps", m_minsubsteps, 1);
60  m_session->LoadParameter("MaxSubSteps", m_maxsubsteps, 100);
61 
62  int dim = m_fields[0]->GetCoordim(0);
64 
65  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
66  for (int i = 0; i < dim; ++i)
67  {
68  m_traceNormals[i] = Array<OneD, NekDouble>(nTracePts);
69  }
70  m_fields[0]->GetTrace()->GetNormals(m_traceNormals);
71 }
72 
74 {
75 }
76 
78  const Array<OneD, const Array<OneD, NekDouble>> &fields,
79  const Array<OneD, const Array<OneD, NekDouble>> &N, NekDouble kinvis)
80 {
81  ASSERTL0(false, "This method should not be called by Substepping routine");
82 }
83 
85  const LibUtilities::TimeIntegrationSchemeSharedPtr &IntegrationScheme)
86 {
87  m_intScheme = IntegrationScheme;
88 
89  unsigned int order = IntegrationScheme->GetOrder();
90 
91  // Set to 1 for first step and it will then be increased in
92  // time advance routines
93  if ((IntegrationScheme->GetName() == "Euler" &&
94  IntegrationScheme->GetVariant() == "Backward") ||
95 
96  (IntegrationScheme->GetName() == "BDFImplicit" &&
97  (order == 1 || order == 2)))
98  {
99  // Note RK first order SSP is just Forward Euler.
100  std::string vSubStepIntScheme = "RungeKutta";
101  std::string vSubStepIntSchemeVariant = "SSP";
102  int vSubStepIntSchemeOrder = order;
103 
104  if (m_session->DefinesSolverInfo("SubStepIntScheme"))
105  {
106  vSubStepIntScheme = m_session->GetSolverInfo("SubStepIntScheme");
107  vSubStepIntSchemeVariant = "";
108  vSubStepIntSchemeOrder = order;
109  }
110 
113  vSubStepIntScheme, vSubStepIntSchemeVariant,
114  vSubStepIntSchemeOrder, std::vector<NekDouble>());
115 
116  int nvel = m_velocity.size();
117  int ndim = order + 1;
118 
119  // Fields for linear/quadratic interpolation
121  int ntotpts = m_fields[0]->GetTotPoints();
122  m_previousVelFields[0] = Array<OneD, NekDouble>(ndim * nvel * ntotpts);
123 
124  for (int i = 1; i < ndim * nvel; ++i)
125  {
126  m_previousVelFields[i] = m_previousVelFields[i - 1] + ntotpts;
127  }
128  }
129  else
130  {
131  ASSERTL0(0, "Integration method not suitable: Options include "
132  "BackwardEuler or BDFImplicitOrder{1,2}");
133  }
134 
135  m_intSteps = IntegrationScheme->GetNumIntegrationPhases();
136 
137  // set explicit time-integration class operators
142 }
143 
144 /**
145  * Explicit Advection terms used by SubStepAdvance time integration
146  */
148  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
149  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
150 {
151  int i;
152  int nVariables = inarray.size();
153  int nQuadraturePts = inarray[0].size();
154 
155  /// Get the number of coefficients
156  int ncoeffs = m_fields[0]->GetNcoeffs();
157 
158  /// Define an auxiliary variable to compute the RHS
159  Array<OneD, Array<OneD, NekDouble>> WeakAdv(nVariables);
160  WeakAdv[0] = Array<OneD, NekDouble>(ncoeffs * nVariables);
161  for (i = 1; i < nVariables; ++i)
162  {
163  WeakAdv[i] = WeakAdv[i - 1] + ncoeffs;
164  }
165 
167 
168  Velfields[0] = Array<OneD, NekDouble>(nQuadraturePts * m_velocity.size());
169 
170  for (i = 1; i < m_velocity.size(); ++i)
171  {
172  Velfields[i] = Velfields[i - 1] + nQuadraturePts;
173  }
174 
175  SubStepExtrapolateField(fmod(time, m_timestep), Velfields);
176 
177  for (auto &x : m_forcing)
178  {
179  x->PreApply(m_fields, Velfields, Velfields, time);
180  }
181  m_advObject->Advect(m_velocity.size(), m_fields, Velfields, inarray,
182  outarray, time);
183  for (auto &x : m_forcing)
184  {
185  x->Apply(m_fields, outarray, outarray, time);
186  }
187 
188  for (i = 0; i < nVariables; ++i)
189  {
190  m_fields[i]->IProductWRTBase(outarray[i], WeakAdv[i]);
191  // negation requried due to sign of DoAdvection term to be consistent
192  Vmath::Neg(ncoeffs, WeakAdv[i], 1);
193  }
194 
195  AddAdvectionPenaltyFlux(Velfields, inarray, WeakAdv);
196 
197  /// Operations to compute the RHS
198  for (i = 0; i < nVariables; ++i)
199  {
200  // Negate the RHS
201  Vmath::Neg(ncoeffs, WeakAdv[i], 1);
202 
203  /// Multiply the flux by the inverse of the mass matrix
204  m_fields[i]->MultiplyByElmtInvMass(WeakAdv[i], WeakAdv[i]);
205 
206  /// Store in outarray the physical values of the RHS
207  m_fields[i]->BwdTrans(WeakAdv[i], outarray[i]);
208  }
209 }
210 
211 /**
212  * Projection used by SubStepAdvance time integration
213  */
215  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
216  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
217 {
218  ASSERTL1(inarray.size() == outarray.size(),
219  "Inarray and outarray of different sizes ");
220 
221  for (int i = 0; i < inarray.size(); ++i)
222  {
223  Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
224  }
225 }
226 
227 /**
228  *
229  */
231  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
232  const NekDouble Aii_Dt, NekDouble kinvis)
233 {
234  // int nConvectiveFields =m_fields.size()-1;
235  Array<OneD, Array<OneD, NekDouble>> nullvelfields;
236 
237  m_pressureCalls++;
238 
239  // Calculate non-linear and viscous BCs at current level and
240  // put in m_pressureHBCs[0]
241  CalcNeumannPressureBCs(inarray, nullvelfields, kinvis);
242 
243  // Extrapolate to m_pressureHBCs to n+1
245 
246  // Add (phi,Du/Dt) term to m_presureHBC
247  AddDuDt();
248 
249  // Copy m_pressureHBCs to m_PbndExp
251 
252  // Evaluate High order outflow conditiosn if required.
253  CalcOutflowBCs(inarray, kinvis);
254 }
255 
256 /**
257  *
258  */
260 {
261  int i, n;
262  int nvel = m_velocity.size();
263  int npts = m_fields[0]->GetTotPoints();
264 
265  // rotate fields
266  int nblocks = m_previousVelFields.size() / nvel;
268 
269  // rotate storage space
270  for (n = 0; n < nvel; ++n)
271  {
272  save = m_previousVelFields[(nblocks - 1) * nvel + n];
273 
274  for (i = nblocks - 1; i > 0; --i)
275  {
276  m_previousVelFields[i * nvel + n] =
277  m_previousVelFields[(i - 1) * nvel + n];
278  }
279 
280  m_previousVelFields[n] = save;
281  }
282 
283  // Put previous field
284  for (i = 0; i < nvel; ++i)
285  {
286  m_fields[m_velocity[i]]->BwdTrans(
287  m_fields[m_velocity[i]]->GetCoeffs(),
288  m_fields[m_velocity[i]]->UpdatePhys());
289  Vmath::Vcopy(npts, m_fields[m_velocity[i]]->GetPhys(), 1,
290  m_previousVelFields[i], 1);
291  }
292 
293  if (nstep == 0) // initialise all levels with first field
294  {
295  for (n = 0; n < nvel; ++n)
296  {
297  for (i = 1; i < nblocks; ++i)
298  {
299  Vmath::Vcopy(npts, m_fields[m_velocity[n]]->GetPhys(), 1,
300  m_previousVelFields[i * nvel + n], 1);
301  }
302  }
303  }
304 }
305 
306 /**
307  *
308  */
310 {
311  int n;
312  int nsubsteps;
313 
314  NekDouble dt;
315 
317 
318  static int ncalls = 1;
319  int nint = std::min(ncalls++, m_intSteps);
320 
321  // this needs to change
322  m_comm = m_fields[0]->GetComm()->GetRowComm();
323 
324  // Get the proper time step with CFL control
325  dt = GetSubstepTimeStep();
326 
327  nsubsteps = (m_timestep > dt) ? ((int)(m_timestep / dt) + 1) : 1;
328  nsubsteps = std::max(m_minsubsteps, nsubsteps);
329 
330  ASSERTL0(nsubsteps < m_maxsubsteps,
331  "Number of substeps has exceeded maximum");
332 
333  dt = m_timestep / nsubsteps;
334 
335  if (m_infosteps && !((nstep + 1) % m_infosteps) && m_comm->GetRank() == 0)
336  {
337  std::cout << "Sub-integrating using " << nsubsteps
338  << " steps over Dt = " << m_timestep
339  << " (SubStep CFL=" << m_cflSafetyFactor << ")" << std::endl;
340  }
341 
342  const TripleArray &solutionVector = m_intScheme->GetSolutionVector();
343 
344  for (int m = 0; m < nint; ++m)
345  {
346  // We need to update the fields held by the m_intScheme
347  fields = solutionVector[m];
348 
349  // Initialise NS solver which is set up to use a GLM method
350  // with calls to EvaluateAdvection_SetPressureBCs and
351  // SolveUnsteadyStokesSystem
352  m_subStepIntegrationScheme->InitializeScheme(dt, fields, time,
354 
355  for (n = 0; n < nsubsteps; ++n)
356  {
357  fields = m_subStepIntegrationScheme->TimeIntegrate(
358  n, dt, m_subStepIntegrationOps);
359  }
360 
361  // set up HBC m_acceleration field for Pressure BCs
363 
364  // Reset time integrated solution in m_intScheme
365  m_intScheme->SetSolutionVector(m, fields);
366  }
367 }
368 
369 /**
370  *
371  */
373 {
374  int n_element = m_fields[0]->GetExpSize();
375 
376  const Array<OneD, int> ExpOrder = m_fields[0]->EvalBasisNumModesMaxPerExp();
377 
378  const NekDouble cLambda = 0.2; // Spencer book pag. 317
379 
380  Array<OneD, NekDouble> tstep(n_element, 0.0);
381  Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
383 
384  for (int i = 0; i < m_velocity.size(); ++i)
385  {
386  velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
387  }
388  stdVelocity = GetMaxStdVelocity(velfields);
389 
390  for (int el = 0; el < n_element; ++el)
391  {
392  tstep[el] =
393  m_cflSafetyFactor / (stdVelocity[el] * cLambda *
394  (ExpOrder[el] - 1) * (ExpOrder[el] - 1));
395  }
396 
397  NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
398  m_comm->AllReduce(TimeStep, LibUtilities::ReduceMin);
399 
400  return TimeStep;
401 }
402 
404  const Array<OneD, const Array<OneD, NekDouble>> &velfield,
405  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
406  Array<OneD, Array<OneD, NekDouble>> &Outarray)
407 {
408  ASSERTL1(physfield.size() == Outarray.size(),
409  "Physfield and outarray are of different dimensions");
410 
411  int i;
412 
413  /// Number of trace points
414  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
415 
416  /// Number of spatial dimensions
417  int nDimensions = m_bnd_dim;
418 
419  /// Forward state array
420  Array<OneD, NekDouble> Fwd(3 * nTracePts);
421 
422  /// Backward state array
423  Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
424 
425  /// upwind numerical flux state array
426  Array<OneD, NekDouble> numflux = Bwd + nTracePts;
427 
428  /// Normal velocity array
429  Array<OneD, NekDouble> Vn(nTracePts, 0.0);
430 
431  // Extract velocity field along the trace space and multiply by trace
432  // normals
433  for (i = 0; i < nDimensions; ++i)
434  {
435  m_fields[0]->ExtractTracePhys(velfield[i], Fwd);
436  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Fwd, 1, Vn, 1, Vn, 1);
437  }
438 
439  for (i = 0; i < physfield.size(); ++i)
440  {
441  /// Extract forwards/backwards trace spaces
442  /// Note: Needs to have correct i value to get boundary conditions
443  m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
444 
445  /// Upwind between elements
446  m_fields[0]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux);
447 
448  /// Construct difference between numflux and Fwd,Bwd
449  Vmath::Vsub(nTracePts, numflux, 1, Fwd, 1, Fwd, 1);
450  Vmath::Vsub(nTracePts, numflux, 1, Bwd, 1, Bwd, 1);
451 
452  /// Calculate the numerical fluxes multipling Fwd, Bwd and
453  /// numflux by the normal advection velocity
454  Vmath::Vmul(nTracePts, Fwd, 1, Vn, 1, Fwd, 1);
455  Vmath::Vmul(nTracePts, Bwd, 1, Vn, 1, Bwd, 1);
456 
457  m_fields[0]->AddFwdBwdTraceIntegral(Fwd, Bwd, Outarray[i]);
458  }
459 }
460 
461 /**
462  * Extrapolate field using equally time spaced field un,un-1,un-2, (at
463  * dt intervals) to time n+t at order Ord
464  */
466  NekDouble toff, Array<OneD, Array<OneD, NekDouble>> &ExtVel)
467 {
468  int npts = m_fields[0]->GetTotPoints();
469  int nvel = m_velocity.size();
470  int i, j;
472 
473  int ord = m_intSteps;
474 
475  // calculate Lagrange interpolants
476  Vmath::Fill(4, 1.0, l, 1);
477 
478  for (i = 0; i <= ord; ++i)
479  {
480  for (j = 0; j <= ord; ++j)
481  {
482  if (i != j)
483  {
484  l[i] *= (j * m_timestep + toff);
485  l[i] /= (j * m_timestep - i * m_timestep);
486  }
487  }
488  }
489 
490  for (i = 0; i < nvel; ++i)
491  {
492  Vmath::Smul(npts, l[0], m_previousVelFields[i], 1, ExtVel[i], 1);
493 
494  for (j = 1; j <= ord; ++j)
495  {
496  Blas::Daxpy(npts, l[j], m_previousVelFields[j * nvel + i], 1,
497  ExtVel[i], 1);
498  }
499  }
500 }
501 
502 /**
503  *
504  */
506  int HBCdata, NekDouble kinvis, Array<OneD, NekDouble> &Q,
508 {
509  Vmath::Smul(HBCdata, -kinvis, Q, 1, Q, 1);
510 }
511 
513 {
514  return m_subStepIntegrationScheme->GetName();
515 }
516 
517 } // 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
Array< OneD, Array< OneD, NekDouble > > m_pressureHBCs
Storage for current and previous levels of high order pressure boundary conditions.
Definition: Extrapolate.h:241
int m_bnd_dim
bounday dimensionality
Definition: Extrapolate.h:217
Array< OneD, NekDouble > GetMaxStdVelocity(const Array< OneD, Array< OneD, NekDouble >> inarray)
void CalcNeumannPressureBCs(const Array< OneD, const Array< OneD, NekDouble >> &fields, const Array< OneD, const Array< OneD, NekDouble >> &N, NekDouble kinvis)
Definition: Extrapolate.h:170
void CopyPressureHBCsToPbndExp(void)
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Velocity fields.
Definition: Extrapolate.h:198
Array< OneD, Array< OneD, NekDouble > > m_iprodnormvel
Storage for current and previous levels of the inner product of normal velocity.
Definition: Extrapolate.h:245
void ExtrapolateArray(Array< OneD, Array< OneD, NekDouble >> &array)
int m_intSteps
Maximum points used in pressure BC evaluation.
Definition: Extrapolate.h:235
void CalcOutflowBCs(const Array< OneD, const Array< OneD, NekDouble >> &fields, NekDouble kinvis)
NekDouble m_timestep
Definition: Extrapolate.h:237
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Definition: Extrapolate.h:209
SolverUtils::AdvectionSharedPtr m_advObject
Definition: Extrapolate.h:207
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
Definition: Extrapolate.h:205
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: Extrapolate.h:247
int m_pressureCalls
number of times the high-order pressure BCs have been called
Definition: Extrapolate.h:226
LibUtilities::SessionReaderSharedPtr m_session
Definition: Extrapolate.h:190
void IProductNormVelocityOnHBC(const Array< OneD, const Array< OneD, NekDouble >> &Vel, Array< OneD, NekDouble > &IprodVn)
LibUtilities::CommSharedPtr m_comm
Definition: Extrapolate.h:192
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)
An abstract base class encapsulating the concept of advection of a vector field.
Definition: Advection.h:72
void SubStepExtrapolateField(NekDouble toff, Array< OneD, Array< OneD, NekDouble >> &ExtVel)
virtual void v_SubStepAdvance(int nstep, NekDouble time)
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps
virtual std::string v_GetSubStepName(void)
LibUtilities::TimeIntegrationSchemeSharedPtr m_subStepIntegrationScheme
static std::string className
Name of class.
static ExtrapolateSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, MultiRegions::ExpListSharedPtr &pPressure, const Array< OneD, int > &pVel, const SolverUtils::AdvectionSharedPtr &advObject)
Creates an instance of this class.
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
void AddAdvectionPenaltyFlux(const Array< OneD, const Array< OneD, NekDouble >> &velfield, const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &outarray)
void SubStepAdvection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
virtual void v_EvaluatePressureBCs(const Array< OneD, const Array< OneD, NekDouble >> &fields, const Array< OneD, const Array< OneD, NekDouble >> &N, NekDouble kinvis)
virtual void v_SubSteppingTimeIntegration(const LibUtilities::TimeIntegrationSchemeSharedPtr &IntegrationScheme)
virtual void v_SubStepSaveFields(int nstep)
SubSteppingExtrapolate(const LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields, MultiRegions::ExpListSharedPtr pPressure, const Array< OneD, int > pVel, const SolverUtils::AdvectionSharedPtr advObject)
virtual void v_MountHOPBCs(int HBCdata, NekDouble kinvis, Array< OneD, NekDouble > &Q, Array< OneD, const NekDouble > &Advection)
void SubStepProjection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
virtual void v_SubStepSetPressureBCs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, NekDouble Aii_Dt, NekDouble kinvis)
Array< OneD, Array< OneD, NekDouble > > m_previousVelFields
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:154
TimeIntegrationSchemeFactory & GetTimeIntegrationSchemeFactory()
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:280
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
ExtrapolateFactory & GetExtrapolateFactory()
Definition: Extrapolate.cpp:48
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
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.cpp:1050
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 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 Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:419