Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
VCSMapping.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File VCSMapping.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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Velocity Correction Scheme w/ coordinate transformation
33 // for the Incompressible Navier Stokes equations
34 ///////////////////////////////////////////////////////////////////////////////
35 
38 #include <SolverUtils/Core/Misc.h>
39 
40 #include <boost/algorithm/string.hpp>
41 
42 namespace Nektar
43 {
44  string VCSMapping::className =
46  "VCSMapping",
48 
49  /**
50  * Constructor. Creates ...
51  *
52  * \param
53  * \param
54  */
57  : UnsteadySystem(pSession),
58  VelocityCorrectionScheme(pSession)
59  {
60 
61  }
62 
64  {
66 
69  "Could not create mapping in VCSMapping.");
70 
71  std::string vExtrapolation = "Mapping";
73  vExtrapolation,
74  m_session,
75  m_fields,
76  m_pressure,
77  m_velocity,
78  m_advObject);
79  m_extrapolation->SubSteppingTimeIntegration(
80  m_intScheme->GetIntegrationMethod(), m_intScheme);
81  m_extrapolation->GenerateHOPBCMap();
82 
83  // Storage to extrapolate pressure forcing
84  int physTot = m_fields[0]->GetTotPoints();
85  int intSteps = 1;
86  int intMethod = m_intScheme->GetIntegrationMethod();
87  switch(intMethod)
88  {
90  {
91  intSteps = 1;
92  }
93  break;
95  {
96  intSteps = 2;
97  }
98  break;
100  {
101  intSteps = 3;
102  }
103  break;
104  }
106  for(int i = 0; i < m_presForcingCorrection.num_elements(); i++)
107  {
109  }
110  m_verbose = (m_session->DefinesCmdLineArgument("verbose"))? true :false;
111 
112  // Load solve parameters related to the mapping
113  // Flags determining if pressure/viscous terms should be treated implicitly
114  m_session->MatchSolverInfo("MappingImplicitPressure","True",
115  m_implicitPressure,false);
116  m_session->MatchSolverInfo("MappingImplicitViscous","True",
117  m_implicitViscous,false);
118  m_session->MatchSolverInfo("MappingNeglectViscous","True",
119  m_neglectViscous,false);
120 
121  if (m_neglectViscous)
122  {
123  m_implicitViscous = false;
124  }
125 
126  // Tolerances and relaxation parameters for implicit terms
127  m_session->LoadParameter("MappingPressureTolerance",
128  m_pressureTolerance,1e-12);
129  m_session->LoadParameter("MappingViscousTolerance",
130  m_viscousTolerance,1e-12);
131  m_session->LoadParameter("MappingPressureRelaxation",
133  m_session->LoadParameter("MappingViscousRelaxation",
134  m_viscousRelaxation,1.0);
135 
136  }
137 
138  /**
139  * Destructor
140  */
142  {
143  }
144 
146  {
147  UnsteadySystem::v_DoInitialise();
148 
149  // Set up Field Meta Data for output files
150  m_fieldMetaDataMap["Kinvis"] =
151  boost::lexical_cast<std::string>(m_kinvis);
152  m_fieldMetaDataMap["TimeStep"] =
153  boost::lexical_cast<std::string>(m_timestep);
154 
155  // Correct Dirichlet boundary conditions to account for mapping
156  m_mapping->UpdateBCs(0.0);
157  //
158  for(int i = 0; i < m_nConvectiveFields; ++i)
159  {
160  m_fields[i]->LocalToGlobal();
161  m_fields[i]->ImposeDirichletConditions(m_fields[i]->UpdateCoeffs());
162  m_fields[i]->GlobalToLocal();
163  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
164  m_fields[i]->UpdatePhys());
165  }
166 
167  // Initialise m_gradP
168  int physTot = m_fields[0]->GetTotPoints();
170  for(int i = 0; i < m_nConvectiveFields; ++i)
171  {
172  m_gradP[i] = Array<OneD, NekDouble>(physTot,0.0);
174  m_pressure->GetPhys(),
175  m_gradP[i]);
176  if(m_pressure->GetWaveSpace())
177  {
178  m_pressure->HomogeneousBwdTrans(m_gradP[i],m_gradP[i]);
179  }
180  }
181  }
182 
183  /**
184  * Explicit part of the method - Advection, Forcing + HOPBCs
185  */
187  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
188  Array<OneD, Array<OneD, NekDouble> > &outarray,
189  const NekDouble time)
190  {
191  // Update mapping and Deal with Dirichlet boundary conditions
192  if (m_mapping->IsTimeDependent())
193  {
194  if (m_mapping->IsFromFunction())
195  {
196  // If the transformation is explicitly defined, update it here
197  // Otherwise, it will be done somewhere else (ForcingMovingBody)
198  m_mapping->UpdateMapping(time);
199  }
200  m_mapping->UpdateBCs(time);
201  }
202 
203  EvaluateAdvectionTerms(inarray, outarray);
204 
205  // Smooth advection
207  {
208  for(int i = 0; i < m_nConvectiveFields; ++i)
209  {
210  m_pressure->SmoothField(outarray[i]);
211  }
212  }
213 
214  // Add forcing terms
215  std::vector<SolverUtils::ForcingSharedPtr>::const_iterator x;
216  for (x = m_forcing.begin(); x != m_forcing.end(); ++x)
217  {
218  (*x)->Apply(m_fields, inarray, outarray, time);
219  }
220 
221  // Add mapping terms
222  ApplyIncNSMappingForcing( outarray );
223 
224  // Calculate High-Order pressure boundary conditions
225  m_extrapolation->EvaluatePressureBCs(inarray,outarray,m_kinvis);
226  }
227 
228 
229  /**
230  * Forcing term for Poisson solver solver
231  */
233  const Array<OneD, const Array<OneD, NekDouble> > &fields,
234  Array<OneD, Array<OneD, NekDouble> > &Forcing,
235  const NekDouble aii_Dt)
236  {
237  if (m_mapping->HasConstantJacobian())
238  {
240  Forcing, aii_Dt);
241  }
242  else
243  {
244  int physTot = m_fields[0]->GetTotPoints();
245  int nvel = m_nConvectiveFields;
246  Array<OneD, NekDouble> wk(physTot, 0.0);
247 
248  Array<OneD, NekDouble> Jac(physTot,0.0);
249  m_mapping->GetJacobian(Jac);
250 
251  // Calculate div(J*u/Dt)
252  Vmath::Zero(physTot,Forcing[0],1);
253  for(int i = 0; i < nvel; ++i)
254  {
255  if (m_fields[i]->GetWaveSpace())
256  {
257  m_fields[i]->HomogeneousBwdTrans(fields[i],wk);
258  }
259  else
260  {
261  Vmath::Vcopy(physTot, fields[i], 1, wk, 1);
262  }
263  Vmath::Vmul(physTot,wk,1,Jac,1,wk,1);
264  if (m_fields[i]->GetWaveSpace())
265  {
266  m_fields[i]->HomogeneousFwdTrans(wk,wk);
267  }
268  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[i],wk, wk);
269  Vmath::Vadd(physTot,wk,1,Forcing[0],1,Forcing[0],1);
270  }
271  Vmath::Smul(physTot,1.0/aii_Dt,Forcing[0],1,Forcing[0],1);
272 
273  //
274  // If the mapping viscous terms are being treated explicitly
275  // we need to apply a correction to the forcing
276  if (!m_implicitViscous)
277  {
278  bool wavespace = m_fields[0]->GetWaveSpace();
279  m_fields[0]->SetWaveSpace(false);
280 
281  //
282  // Part 1: div(J*grad(U/J . grad(J)))
284  Array<OneD, Array<OneD, NekDouble> > velocity (nvel);
285  for(int i = 0; i < tmp.num_elements(); i++)
286  {
287  tmp[i] = Array<OneD, NekDouble>(physTot,0.0);
288  velocity[i] = Array<OneD, NekDouble>(physTot,0.0);
289  if (wavespace)
290  {
291  m_fields[0]->HomogeneousBwdTrans(m_fields[i]->GetPhys(),
292  velocity[i]);
293  }
294  else
295  {
296  Vmath::Vcopy(physTot, m_fields[i]->GetPhys(), 1,
297  velocity[i], 1);
298  }
299  }
300  // Calculate wk = U.grad(J)
301  m_mapping->DotGradJacobian(velocity, wk);
302  // Calculate wk = (U.grad(J))/J
303  Vmath::Vdiv(physTot, wk, 1, Jac, 1, wk, 1);
304  // J*grad[(U.grad(J))/J]
305  for(int i = 0; i < nvel; ++i)
306  {
307  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i],
308  wk, tmp[i]);
309  Vmath::Vmul(physTot, Jac, 1, tmp[i], 1, tmp[i], 1);
310  }
311  // div(J*grad[(U.grad(J))/J])
312  Vmath::Zero(physTot, wk, 1);
313  for(int i = 0; i < nvel; ++i)
314  {
315  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i],
316  tmp[i], tmp[i]);
317  Vmath::Vadd(physTot, wk, 1, tmp[i], 1, wk, 1);
318  }
319 
320  // Part 2: grad(J) . curl(curl(U))
321  m_mapping->CurlCurlField(velocity, tmp, m_implicitViscous);
322  // dont need velocity any more, so reuse it
323  m_mapping->DotGradJacobian(tmp, velocity[0]);
324 
325  // Add two parts
326  Vmath::Vadd(physTot, velocity[0], 1, wk, 1, wk, 1);
327 
328  // Multiply by kinvis
329  Vmath::Smul(physTot, m_kinvis, wk, 1, wk, 1);
330 
331  // Extrapolate correction
332  m_extrapolation->ExtrapolateArray(m_presForcingCorrection,
333  wk, wk);
334 
335  // Put in wavespace
336  if (wavespace)
337  {
338  m_fields[0]->HomogeneousFwdTrans(wk,wk);
339  }
340  // Apply correction: Forcing = Forcing - correction
341  Vmath::Vsub(physTot, Forcing[0], 1, wk, 1, Forcing[0], 1);
342 
343  m_fields[0]->SetWaveSpace(wavespace);
344  }
345  }
346  }
347 
348  /**
349  * Forcing term for Helmholtz solver
350  */
352  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
353  Array<OneD, Array<OneD, NekDouble> > &Forcing,
354  const NekDouble aii_Dt)
355  {
356  NekDouble aii_dtinv = 1.0/aii_Dt;
357  int physTot = m_fields[0]->GetTotPoints();
358 
359  // Grad p
360  m_pressure->BwdTrans(m_pressure->GetCoeffs(),m_pressure->UpdatePhys());
361 
362  int nvel = m_velocity.num_elements();
363  if(nvel == 2)
364  {
365  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1]);
366  }
367  else
368  {
369  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1],
370  Forcing[2]);
371  }
372 
373  // Copy grad p in physical space to m_gradP to reuse later
374  if (m_pressure->GetWaveSpace())
375  {
376  for (int i=0; i<nvel; i++)
377  {
378  m_pressure->HomogeneousBwdTrans(Forcing[i],m_gradP[i]);
379  }
380  }
381  else
382  {
383  for (int i=0; i<nvel; i++)
384  {
385  Vmath::Vcopy(physTot, Forcing[i], 1, m_gradP[i], 1);
386  }
387  }
388 
389  if ( (!m_mapping->HasConstantJacobian()) || m_implicitPressure)
390  {
391  // If pressure terms are treated explicitly, we need to divide by J
392  // if they are implicit, we need to calculate G(p)
393  if (m_implicitPressure)
394  {
395  m_mapping->RaiseIndex(m_gradP, Forcing);
396  }
397  else
398  {
399  Array<OneD, NekDouble> Jac(physTot,0.0);
400  m_mapping->GetJacobian(Jac);
401  for (int i=0; i<nvel; i++)
402  {
403  Vmath::Vdiv(physTot, m_gradP[i], 1, Jac, 1, Forcing[i], 1);
404  }
405  }
406  // Transform back to wavespace
407  if (m_pressure->GetWaveSpace())
408  {
409  for (int i=0; i<nvel; i++)
410  {
411  m_pressure->HomogeneousFwdTrans(Forcing[i],Forcing[i]);
412  }
413  }
414  }
415 
416  // Subtract inarray/(aii_dt) and divide by kinvis. Kinvis will
417  // need to be updated for the convected fields.
418  for(int i = 0; i < nvel; ++i)
419  {
420  Blas::Daxpy(physTot,-aii_dtinv,inarray[i],1,Forcing[i],1);
421  Blas::Dscal(physTot,1.0/m_kinvis,&(Forcing[i])[0],1);
422  }
423  }
424 
425  /**
426  * Solve pressure system
427  */
429  const Array<OneD, NekDouble> &Forcing)
430  {
431  if (!m_implicitPressure)
432  {
434  }
435  else
436  {
437  int physTot = m_fields[0]->GetTotPoints();
438  int nvel = m_nConvectiveFields;
439  bool converged = false; // flag to mark if system converged
440  int s = 0; // iteration counter
441  NekDouble error; // L2 error at current iteration
442  NekDouble forcing_L2 = 0.0; // L2 norm of F
443 
444  int maxIter;
445  m_session->LoadParameter("MappingMaxIter",maxIter,5000);
446 
447  // rhs of the equation at current iteration
448  Array< OneD, NekDouble> F_corrected(physTot, 0.0);
449  // Pressure field at previous iteration
450  Array<OneD, NekDouble> previous_iter (physTot, 0.0);
451  // Temporary variables
455  for(int i = 0; i < nvel; ++i)
456  {
457  wk1[i] = Array<OneD, NekDouble> (physTot, 0.0);
458  wk2[i] = Array<OneD, NekDouble> (physTot, 0.0);
459  gradP[i] = Array<OneD, NekDouble> (physTot, 0.0);
460  }
461 
462  // Jacobian
463  Array<OneD, NekDouble> Jac(physTot, 0.0);
464  m_mapping->GetJacobian(Jac);
465 
466  // Factors for Laplacian system
468  factors[StdRegions::eFactorLambda] = 0.0;
469 
470  m_pressure->BwdTrans(m_pressure->GetCoeffs(),
471  m_pressure->UpdatePhys());
472  forcing_L2 = m_pressure->L2(Forcing, wk1[0]);
473  while (!converged)
474  {
475  // Update iteration counter and set previous iteration field
476  // (use previous timestep solution for first iteration)
477  s++;
478  ASSERTL0(s < maxIter,
479  "VCSMapping exceeded maximum number of iterations.");
480 
481  Vmath::Vcopy(physTot, m_pressure->GetPhys(), 1,
482  previous_iter, 1);
483 
484  // Correct pressure bc to account for iteration
485  m_extrapolation->CorrectPressureBCs(previous_iter);
486 
487  //
488  // Calculate forcing term for this iteration
489  //
490  for(int i = 0; i < nvel; ++i)
491  {
493  previous_iter, gradP[i]);
494  if(m_pressure->GetWaveSpace())
495  {
496  m_pressure->HomogeneousBwdTrans(gradP[i], wk1[i]);
497  }
498  else
499  {
500  Vmath::Vcopy(physTot, gradP[i], 1, wk1[i], 1);
501  }
502  }
503  m_mapping->RaiseIndex(wk1, wk2); // G(p)
504 
505  m_mapping->Divergence(wk2, F_corrected); // div(G(p))
506  if (!m_mapping->HasConstantJacobian())
507  {
508  Vmath::Vmul(physTot, F_corrected, 1,
509  Jac, 1,
510  F_corrected, 1);
511  }
512  // alpha*J*div(G(p))
513  Vmath::Smul(physTot, m_pressureRelaxation, F_corrected, 1,
514  F_corrected, 1);
515  if(m_pressure->GetWaveSpace())
516  {
517  m_pressure->HomogeneousFwdTrans(F_corrected, F_corrected);
518  }
519  // alpha*J*div(G(p)) - p_ii
520  for (int i = 0; i < m_nConvectiveFields; ++i)
521  {
523  gradP[i], wk1[0]);
524  Vmath::Vsub(physTot, F_corrected, 1, wk1[0], 1,
525  F_corrected, 1);
526  }
527  // p_i,i - J*div(G(p))
528  Vmath::Neg(physTot, F_corrected, 1);
529  // alpha*F - alpha*J*div(G(p)) + p_i,i
530  Vmath::Smul(physTot, m_pressureRelaxation, Forcing, 1,
531  wk1[0], 1);
532  Vmath::Vadd(physTot, wk1[0], 1, F_corrected, 1, F_corrected, 1);
533 
534  //
535  // Solve system
536  //
537  m_pressure->HelmSolve(F_corrected, m_pressure->UpdateCoeffs(),
538  NullFlagList, factors);
539  m_pressure->BwdTrans(m_pressure->GetCoeffs(),
540  m_pressure->UpdatePhys());
541 
542  //
543  // Test convergence
544  //
545  error = m_pressure->L2(m_pressure->GetPhys(), previous_iter);
546  if ( forcing_L2 != 0)
547  {
548  if ( (error/forcing_L2 < m_pressureTolerance))
549  {
550  converged = true;
551  }
552  }
553  else
554  {
555  if ( error < m_pressureTolerance)
556  {
557  converged = true;
558  }
559  }
560  }
561  if (m_verbose && m_session->GetComm()->GetRank()==0)
562  {
563  std::cout << " Pressure system (mapping) converged in " << s <<
564  " iterations with error = " << error << std::endl;
565  }
566  }
567  }
568 
569  /**
570  * Solve velocity system
571  */
573  const Array<OneD, const Array<OneD, NekDouble> > &Forcing,
574  Array<OneD, Array<OneD, NekDouble> > &outarray,
575  const NekDouble aii_Dt)
576  {
577  if(!m_implicitViscous)
578  {
579  VelocityCorrectionScheme::v_SolveViscous(Forcing, outarray, aii_Dt);
580  }
581  else
582  {
583  int physTot = m_fields[0]->GetTotPoints();
584  int nvel = m_nConvectiveFields;
585  bool converged = false; // flag to mark if system converged
586  int s = 0; // iteration counter
587  NekDouble error, max_error; // L2 error at current iteration
588 
589  int maxIter;
590  m_session->LoadParameter("MappingMaxIter",maxIter,5000);
591 
592  //L2 norm of F
594 
595  // rhs of the equation at current iteration
596  Array<OneD, Array<OneD, NekDouble> > F_corrected(nvel);
597  // Solution at previous iteration
598  Array<OneD, Array<OneD, NekDouble> > previous_iter(nvel);
599  // Working space
601  for(int i = 0; i < nvel; ++i)
602  {
603  F_corrected[i] = Array<OneD, NekDouble> (physTot, 0.0);
604  previous_iter[i] = Array<OneD, NekDouble> (physTot, 0.0);
605  wk[i] = Array<OneD, NekDouble> (physTot, 0.0);
606  }
607 
608  // Factors for Helmholtz system
610  factors[StdRegions::eFactorLambda] =
611  1.0*m_viscousRelaxation/aii_Dt/m_kinvis;
612  if(m_useSpecVanVisc)
613  {
617  }
618 
619  // Calculate L2-norm of F and set initial solution for iteration
620  for(int i = 0; i < nvel; ++i)
621  {
622  forcing_L2[i] = m_fields[0]->L2(Forcing[i],wk[0]);
623  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
624  previous_iter[i]);
625  }
626 
627  while (!converged)
628  {
629  converged = true;
630  // Iteration counter
631  s++;
632  ASSERTL0(s < maxIter,
633  "VCSMapping exceeded maximum number of iterations.");
634 
635  max_error = 0.0;
636 
637  //
638  // Calculate forcing term for next iteration
639  //
640 
641  // Calculate L(U)- in this parts all components might be coupled
642  if(m_fields[0]->GetWaveSpace())
643  {
644  for (int i = 0; i < nvel; ++i)
645  {
646  m_fields[0]->HomogeneousBwdTrans(previous_iter[i],
647  wk[i]);
648  }
649  }
650  else
651  {
652  for (int i = 0; i < nvel; ++i)
653  {
654  Vmath::Vcopy(physTot, previous_iter[i], 1, wk[i], 1);
655  }
656  }
657 
658  // (L(U^i) - 1/alpha*U^i_jj)
659  m_mapping->VelocityLaplacian(wk, F_corrected,
660  1.0/m_viscousRelaxation);
661 
662  if(m_fields[0]->GetWaveSpace())
663  {
664  for (int i = 0; i < nvel; ++i)
665  {
666  m_fields[0]->HomogeneousFwdTrans(F_corrected[i],
667  F_corrected[i]);
668  }
669  }
670  else
671  {
672  for (int i = 0; i < nvel; ++i)
673  {
674  Vmath::Vcopy(physTot, F_corrected[i], 1,
675  F_corrected[i], 1);
676  }
677  }
678 
679  // Loop velocity components
680  for (int i = 0; i < nvel; ++i)
681  {
682  // (-alpha*L(U^i) + U^i_jj)
683  Vmath::Smul(physTot, -1.0*m_viscousRelaxation,
684  F_corrected[i], 1,
685  F_corrected[i], 1);
686  // F_corrected = alpha*F + (-alpha*L(U^i) + U^i_jj)
687  Vmath::Smul(physTot, m_viscousRelaxation, Forcing[i], 1,
688  wk[0], 1);
689  Vmath::Vadd(physTot, wk[0], 1, F_corrected[i], 1,
690  F_corrected[i], 1);
691 
692  //
693  // Solve System
694  //
695  m_fields[i]->HelmSolve(F_corrected[i],
696  m_fields[i]->UpdateCoeffs(),
697  NullFlagList, factors);
698  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),outarray[i]);
699 
700  //
701  // Test convergence
702  //
703  error = m_fields[i]->L2(outarray[i], previous_iter[i]);
704 
705  if ( forcing_L2[i] != 0)
706  {
707  if ( (error/forcing_L2[i] >= m_viscousTolerance))
708  {
709  converged = false;
710  }
711  }
712  else
713  {
714  if ( error >= m_viscousTolerance)
715  {
716  converged = false;
717  }
718  }
719  if (error > max_error)
720  {
721  max_error = error;
722  }
723 
724  // Copy field to previous_iter
725  Vmath::Vcopy(physTot, outarray[i], 1, previous_iter[i], 1);
726  }
727  }
728  if (m_verbose && m_session->GetComm()->GetRank()==0)
729  {
730  std::cout << " Velocity system (mapping) converged in " << s <<
731  " iterations with error = " << max_error << std::endl;
732  }
733  }
734  }
735 
736  /**
737  * Explicit terms of the mapping
738  */
740  Array<OneD, Array<OneD, NekDouble> > &outarray)
741  {
742  int physTot = m_fields[0]->GetTotPoints();
747  for (int i = 0; i < m_nConvectiveFields; ++i)
748  {
749  velPhys[i] = Array<OneD, NekDouble> (physTot, 0.0);
750  Forcing[i] = Array<OneD, NekDouble> (physTot, 0.0);
751  tmp[i] = Array<OneD, NekDouble> (physTot, 0.0);
752  }
753 
754  // Get fields and store velocity in wavespace and physical space
755  if(m_fields[0]->GetWaveSpace())
756  {
757  for (int i = 0; i < m_nConvectiveFields; ++i)
758  {
759  vel[i] = m_fields[i]->GetPhys();
760  m_fields[0]->HomogeneousBwdTrans(vel[i],velPhys[i]);
761  }
762  }
763  else
764  {
765  for (int i = 0; i < m_nConvectiveFields; ++i)
766  {
767  vel[i] = m_fields[i]->GetPhys();
768  Vmath::Vcopy(physTot, m_fields[i]->GetPhys(), 1, velPhys[i], 1);
769  }
770  }
771 
772  //Advection contribution
773  MappingAdvectionCorrection(velPhys, Forcing);
774 
775  // Time-derivative contribution
776  if ( m_mapping->IsTimeDependent() )
777  {
778  MappingAccelerationCorrection(vel, velPhys, tmp);
779  for (int i = 0; i < m_nConvectiveFields; ++i)
780  {
781  Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
782  }
783  }
784 
785  // Pressure contribution
786  if (!m_implicitPressure)
787  {
789  for (int i = 0; i < m_nConvectiveFields; ++i)
790  {
791  Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
792  }
793  }
794  // Viscous contribution
795  if ( (!m_implicitViscous) && (!m_neglectViscous))
796  {
797  MappingViscousCorrection(velPhys, tmp);
798  for (int i = 0; i < m_nConvectiveFields; ++i)
799  {
800  Vmath::Smul(physTot, m_kinvis, tmp[i], 1, tmp[i], 1);
801  Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
802  }
803  }
804 
805  // If necessary, transform to wavespace
806  if(m_fields[0]->GetWaveSpace())
807  {
808  for (int i = 0; i < m_nConvectiveFields; ++i)
809  {
810  m_fields[0]->HomogeneousFwdTrans(Forcing[i],Forcing[i]);
811  }
812  }
813 
814  // Add to outarray
815  for (int i = 0; i < m_nConvectiveFields; ++i)
816  {
817  Vmath::Vadd(physTot, outarray[i], 1, Forcing[i], 1, outarray[i], 1);
818  }
819  }
820 
822  const Array<OneD, Array<OneD, NekDouble> > &velPhys,
823  Array<OneD, Array<OneD, NekDouble> > &outarray)
824  {
825  int physTot = m_fields[0]->GetTotPoints();
826  int nvel = m_nConvectiveFields;
827 
828  Array<OneD, Array<OneD, NekDouble> > wk(nvel*nvel);
829 
830  // Apply Christoffel symbols to obtain {i,kj}vel(k)
831  m_mapping->ApplyChristoffelContravar(velPhys, wk);
832 
833  // Calculate correction -U^j*{i,kj}vel(k)
834  for (int i = 0; i< nvel; i++)
835  {
836  Vmath::Zero(physTot,outarray[i],1);
837  for (int j = 0; j< nvel; j++)
838  {
839  Vmath::Vvtvp(physTot,wk[i*nvel+j],1,velPhys[j],1,
840  outarray[i],1,outarray[i],1);
841  }
842  Vmath::Neg(physTot, outarray[i], 1);
843  }
844  }
845 
847  const Array<OneD, Array<OneD, NekDouble> > &vel,
848  const Array<OneD, Array<OneD, NekDouble> > &velPhys,
849  Array<OneD, Array<OneD, NekDouble> > &outarray)
850  {
851  int physTot = m_fields[0]->GetTotPoints();
852  int nvel = m_nConvectiveFields;
853 
854  Array<OneD, Array<OneD, NekDouble> > wk(nvel*nvel);
856  Array<OneD, Array<OneD, NekDouble> > coordVel(nvel);
857  for (int i = 0; i< nvel; i++)
858  {
859  tmp[i] = Array<OneD, NekDouble> (physTot, 0.0);
860  coordVel[i] = Array<OneD, NekDouble> (physTot, 0.0);
861  }
862  // Get coordinates velocity in transformed system
863  m_mapping->GetCoordVelocity(tmp);
864  m_mapping->ContravarFromCartesian(tmp, coordVel);
865 
866  // Calculate first term: U^j u^i,j = U^j (du^i/dx^j + {i,kj}u^k)
867  m_mapping->ApplyChristoffelContravar(velPhys, wk);
868  for (int i=0; i< nvel; i++)
869  {
870  Vmath::Zero(physTot,outarray[i],1);
871 
872  m_fields[0]->PhysDeriv(velPhys[i], tmp[0], tmp[1]);
873  for (int j=0; j< nvel; j++)
874  {
875  if (j == 2)
876  {
877  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j],
878  vel[i], tmp[2]);
879  if (m_fields[0]->GetWaveSpace())
880  {
881  m_fields[0]->HomogeneousBwdTrans(tmp[2],tmp[2]);
882  }
883  }
884 
885  Vmath::Vadd(physTot,wk[i*nvel+j],1,tmp[j],1,
886  wk[i*nvel+j], 1);
887 
888  Vmath::Vvtvp(physTot, coordVel[j], 1, wk[i*nvel+j], 1,
889  outarray[i], 1, outarray[i], 1);
890  }
891  }
892 
893  // Set wavespace to false and store current value
894  bool wavespace = m_fields[0]->GetWaveSpace();
895  m_fields[0]->SetWaveSpace(false);
896 
897  // Add -u^j U^i,j
898  m_mapping->ApplyChristoffelContravar(coordVel, wk);
899  for (int i=0; i< nvel; i++)
900  {
901  if(nvel == 2)
902  {
903  m_fields[0]->PhysDeriv(coordVel[i], tmp[0], tmp[1]);
904  }
905  else
906  {
907  m_fields[0]->PhysDeriv(coordVel[i], tmp[0], tmp[1], tmp[2]);
908  }
909 
910  for (int j=0; j< nvel; j++)
911  {
912  Vmath::Vadd(physTot,wk[i*nvel+j],1,tmp[j],1,
913  wk[i*nvel+j], 1);
914  Vmath::Neg(physTot, wk[i*nvel+j], 1);
915 
916  Vmath::Vvtvp(physTot, velPhys[j], 1, wk[i*nvel+j], 1,
917  outarray[i], 1, outarray[i], 1);
918  }
919  }
920 
921  // Restore value of wavespace
922  m_fields[0]->SetWaveSpace(wavespace);
923  }
924 
926  Array<OneD, Array<OneD, NekDouble> > &outarray)
927  {
928  int physTot = m_fields[0]->GetTotPoints();
929  int nvel = m_nConvectiveFields;
930 
931  // Calculate g^(ij)p_(,j)
932  m_mapping->RaiseIndex(m_gradP, outarray);
933 
934  // Calculate correction = (nabla p)/J - g^(ij)p_,j
935  // (Jac is not required if it is constant)
936  if ( !m_mapping->HasConstantJacobian())
937  {
938  Array<OneD, NekDouble> Jac(physTot, 0.0);
939  m_mapping->GetJacobian(Jac);
940  for(int i = 0; i < nvel; ++i)
941  {
942  Vmath::Vdiv(physTot, m_gradP[i], 1, Jac, 1, m_gradP[i], 1);
943  }
944  }
945  for(int i = 0; i < nvel; ++i)
946  {
947  Vmath::Vsub(physTot, m_gradP[i], 1,outarray[i], 1,
948  outarray[i],1);
949  }
950  }
951 
953  const Array<OneD, Array<OneD, NekDouble> > &velPhys,
954  Array<OneD, Array<OneD, NekDouble> > &outarray)
955  {
956  // L(U) - 1.0*d^2(u^i)/dx^jdx^j
957  m_mapping->VelocityLaplacian(velPhys, outarray, 1.0);
958  }
959 
960 } //end of namespace
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
virtual void v_InitObject()
Init object for UnsteadySystem class.
Definition: VCSMapping.cpp:63
NekDouble m_pressureTolerance
Definition: VCSMapping.h:86
virtual void v_DoInitialise(void)
Sets up initial conditions.
Definition: VCSMapping.cpp:145
void ApplyIncNSMappingForcing(Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:739
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
NekDouble m_kinvis
Kinematic viscosity.
VCSMapping(const LibUtilities::SessionReaderSharedPtr &pSession)
Constructor.
Definition: VCSMapping.cpp:55
ExtrapolateFactory & GetExtrapolateFactory()
Definition: Extrapolate.cpp:48
NekDouble m_timestep
Time step size.
virtual void v_SolveViscous(const Array< OneD, const Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
Definition: VCSMapping.cpp:572
void MappingAccelerationCorrection(const Array< OneD, Array< OneD, NekDouble > > &vel, const Array< OneD, Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:846
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
NekDouble m_sVVDiffCoeff
Diffusion coefficient of SVV modes.
ExtrapolateSharedPtr m_extrapolation
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:75
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:428
virtual void v_SetUpViscousForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
Definition: VCSMapping.cpp:351
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
void Vdiv(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:227
int m_nConvectiveFields
Number of fields to be convected;.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
NekDouble m_viscousRelaxation
Definition: VCSMapping.h:89
void MappingViscousCorrection(const Array< OneD, Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:952
Array< OneD, Array< OneD, NekDouble > > m_gradP
Definition: VCSMapping.h:92
IMEX 2nd order scheme using Backward Different Formula & Extrapolation.
virtual void v_SolvePressure(const Array< OneD, NekDouble > &Forcing)
Definition: VCSMapping.cpp:428
virtual void v_SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
Definition: VCSMapping.cpp:232
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:199
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
bool m_SmoothAdvection
bool to identify if advection term smoothing is requested
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Creates an instance of this class.
Definition: VCSMapping.h:49
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
double NekDouble
virtual void v_SolvePressure(const Array< OneD, NekDouble > &Forcing)
Array< OneD, Array< OneD, NekDouble > > m_presForcingCorrection
Definition: VCSMapping.h:121
virtual ~VCSMapping()
Definition: VCSMapping.cpp:141
void MappingAdvectionCorrection(const Array< OneD, Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:821
EquationSystemFactory & GetEquationSystemFactory()
static std::string className
Name of class.
Definition: VCSMapping.h:59
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:329
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
virtual void v_SolveViscous(const Array< OneD, const Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
void MappingPressureCorrection(Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:925
void EvaluateAdvectionTerms(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
NekDouble m_sVVCutoffRatio
cutt off ratio from which to start decayhing modes
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
bool m_useSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
IMEX 3rd order scheme using Backward Different Formula & Extrapolation.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
NekDouble m_viscousTolerance
Definition: VCSMapping.h:87
virtual void v_InitObject()
Init object for UnsteadySystem class.
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
IMEX 1st order scheme using Euler Backwards/Euler Forwards.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
virtual void v_SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
static GLOBAL_MAPPING_EXPORT MappingSharedPtr Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Return a pointer to the mapping, creating it on first call.
Definition: Mapping.cpp:264
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:285
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:169
virtual void v_EvaluateAdvection_SetPressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: VCSMapping.cpp:186
NekDouble m_pressureRelaxation
Definition: VCSMapping.h:88
static FlagList NullFlagList
An empty flag list.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215