Nektar++
UnsteadySystem.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File UnsteadySystem.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: Generic timestepping for Unsteady solvers
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <iostream>
37 #include <iomanip>
38 
43 
44 namespace Nektar
45 {
46  namespace SolverUtils
47  {
48  /**
49  * @class UnsteadySystem
50  *
51  * Provides the underlying timestepping framework for unsteady solvers
52  * including the general timestepping routines. This class is not
53  * intended to be directly instantiated, but rather is a base class
54  * on which to define unsteady solvers.
55  *
56  * For details on implementing unsteady solvers see
57  * \ref sectionADRSolverModuleImplementation here
58  */
59 
60  /**
61  * Processes SolverInfo parameters from the session file and sets up
62  * timestepping-specific code.
63  * @param pSession Session object to read parameters from.
64  */
67  : EquationSystem(pSession),
68  m_infosteps(10)
69 
70  {
71  }
72 
73  /**
74  * Initialization object for UnsteadySystem class.
75  */
77  {
79 
80  // Load SolverInfo parameters
81  m_session->MatchSolverInfo("DIFFUSIONADVANCEMENT","Explicit",
82  m_explicitDiffusion,true);
83  m_session->MatchSolverInfo("ADVECTIONADVANCEMENT","Explicit",
84  m_explicitAdvection,true);
85  m_session->MatchSolverInfo("REACTIONADVANCEMENT", "Explicit",
86  m_explicitReaction, true);
87 
88  // For steady problems, we do not initialise the time integration
89  if (m_session->DefinesSolverInfo("TIMEINTEGRATIONMETHOD"))
90  {
92  CreateInstance(m_session->GetSolverInfo(
93  "TIMEINTEGRATIONMETHOD"));
94 
95  // Load generic input parameters
96  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
97  m_session->LoadParameter("CFL", m_cflSafetyFactor, 0.0);
98 
99  // Set up time to be dumped in field information
100  m_fieldMetaDataMap["Time"] =
101  boost::lexical_cast<std::string>(m_time);
102  }
103 
104  // By default attempt to forward transform initial condition.
105  m_homoInitialFwd = true;
106 
107  // Set up filters
108  LibUtilities::FilterMap::const_iterator x;
109  LibUtilities::FilterMap f = m_session->GetFilters();
110  for (x = f.begin(); x != f.end(); ++x)
111  {
112  m_filters.push_back(GetFilterFactory().CreateInstance(
113  x->first,
114  m_session,
115  x->second));
116  }
117  }
118 
119  /**
120  * Destructor for the class UnsteadyAdvection.
121  */
123  {
124  }
125 
126  /**
127  * @brief Returns the maximum time estimator for CFL control.
128  */
130  {
131  NekDouble TimeStability = 0.0;
132  switch(m_intScheme->GetIntegrationMethod())
133  {
137  {
138  TimeStability = 2.784;
139  break;
140  }
147  {
148  TimeStability = 2.0;
149  break;
150  }
152  {
153  TimeStability = 1.0;
154  break;
155  }
156  default:
157  {
158  ASSERTL0(
159  false,
160  "No CFL control implementation for this time"
161  "integration scheme");
162  }
163  }
164  return TimeStability;
165  }
166 
167  /**
168  * @brief Initialises the time integration scheme (as specified in the
169  * session file), and perform the time integration.
170  */
172  {
173  ASSERTL0(m_intScheme != 0, "No time integration scheme.");
174 
175  int i, nchk = 1;
176  int nvariables = 0;
177  int nfields = m_fields.num_elements();
178 
179  if (m_intVariables.empty())
180  {
181  for (i = 0; i < nfields; ++i)
182  {
183  m_intVariables.push_back(i);
184  }
185  nvariables = nfields;
186  }
187  else
188  {
189  nvariables = m_intVariables.size();
190  }
191 
192  // Integrate in wave-space if using homogeneous1D
194  {
195  for(i = 0; i < nfields; ++i)
196  {
197  m_fields[i]->HomogeneousFwdTrans(m_fields[i]->GetPhys(),
198  m_fields[i]->UpdatePhys());
199  m_fields[i]->SetWaveSpace(true);
200  m_fields[i]->SetPhysState(false);
201  }
202  }
203 
204  // Set up wrapper to fields data storage.
205  Array<OneD, Array<OneD, NekDouble> > fields(nvariables);
206  Array<OneD, Array<OneD, NekDouble> > tmp (nvariables);
207 
208  // Order storage to list time-integrated fields first.
209  for(i = 0; i < nvariables; ++i)
210  {
211  fields[i] = m_fields[m_intVariables[i]]->GetPhys();
212  m_fields[m_intVariables[i]]->SetPhysState(false);
213  }
214 
215  // Initialise time integration scheme
216  m_intSoln = m_intScheme->InitializeScheme(
217  m_timestep, fields, m_time, m_ode);
218 
219  // Initialise filters
221  for (x = m_filters.begin(); x != m_filters.end(); ++x)
222  {
223  (*x)->Initialise(m_fields, m_time);
224  }
225 
226  // Ensure that there is no conflict of parameters
227  if(m_cflSafetyFactor > 0.0)
228  {
229  // Check final condition
230  ASSERTL0(m_fintime == 0.0 || m_steps == 0,
231  "Final condition not unique: "
232  "fintime > 0.0 and Nsteps > 0");
233 
234  // Check timestep condition
235  ASSERTL0(m_timestep == 0.0,
236  "Timestep not unique: timestep > 0.0 & CFL > 0.0");
237  }
238 
239  // Check uniqueness of checkpoint output
240  ASSERTL0((m_checktime == 0.0 && m_checksteps == 0) ||
241  (m_checktime > 0.0 && m_checksteps == 0) ||
242  (m_checktime == 0.0 && m_checksteps > 0),
243  "Only one of IO_CheckTime and IO_CheckSteps "
244  "should be set!");
245 
246  Timer timer;
247  bool doCheckTime = false;
248  int step = 0;
249  NekDouble intTime = 0.0;
250  NekDouble lastCheckTime = 0.0;
251  NekDouble cpuTime = 0.0;
252  NekDouble elapsed = 0.0;
253 
254  while (step < m_steps ||
256  {
257  if (m_cflSafetyFactor)
258  {
259  m_timestep = GetTimeStep(fields);
260 
261  // Ensure that the final timestep finishes at the final
262  // time, or at a prescribed IO_CheckTime.
263  if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
264  {
266  }
267  else if (m_checktime &&
268  m_time + m_timestep - lastCheckTime >= m_checktime)
269  {
270  lastCheckTime += m_checktime;
271  m_timestep = lastCheckTime - m_time;
272  doCheckTime = true;
273  }
274  }
275 
276  // Perform any solver-specific pre-integration steps
277  if (v_PreIntegrate(step))
278  {
279  break;
280  }
281 
282  timer.Start();
283  fields = m_intScheme->TimeIntegrate(
284  step, m_timestep, m_intSoln, m_ode);
285  timer.Stop();
286 
287  m_time += m_timestep;
288  elapsed = timer.TimePerTest(1);
289  intTime += elapsed;
290  cpuTime += elapsed;
291 
292  // Write out status information
293  if (m_session->GetComm()->GetRank() == 0 &&
294  !((step+1) % m_infosteps))
295  {
296  cout << "Steps: " << setw(8) << left << step+1 << " "
297  << "Time: " << setw(12) << left << m_time;
298 
299  if (m_cflSafetyFactor)
300  {
301  cout << " Time-step: " << setw(12)
302  << left << m_timestep;
303  }
304 
305  stringstream ss;
306  ss << cpuTime << "s";
307  cout << " CPU Time: " << setw(8) << left
308  << ss.str() << endl;
309  cpuTime = 0.0;
310  }
311 
312  // Transform data into coefficient space
313  for (i = 0; i < nvariables; ++i)
314  {
315  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
316  m_fields[m_intVariables[i]]->FwdTrans_IterPerExp(
317  fields[i],
318  m_fields[m_intVariables[i]]->UpdateCoeffs());
319  m_fields[m_intVariables[i]]->SetPhysState(false);
320  }
321 
322  // Perform any solver-specific post-integration steps
323  if (v_PostIntegrate(step))
324  {
325  break;
326  }
327 
328  // search for NaN and quit if found
329  bool nanFound = false;
330  for (i = 0; i < nvariables; ++i)
331  {
332  if (Vmath::Nnan(fields[i].num_elements(), fields[i], 1) > 0)
333  {
334  cout << "NaN found in variable \""
335  << m_session->GetVariable(i)
336  << "\", terminating" << endl;
337  nanFound = true;
338  }
339  }
340 
341  if (nanFound)
342  {
343  break;
344  }
345 
346  // Update filters
348  for (x = m_filters.begin(); x != m_filters.end(); ++x)
349  {
350  (*x)->Update(m_fields, m_time);
351  }
352 
353  // Write out checkpoint files
354  if ((m_checksteps && step && !((step + 1) % m_checksteps)) ||
355  doCheckTime)
356  {
358  {
359  vector<bool> transformed(nfields, false);
360  for(i = 0; i < nfields; i++)
361  {
362  if (m_fields[i]->GetWaveSpace())
363  {
364  m_fields[i]->SetWaveSpace(false);
365  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
366  m_fields[i]->UpdatePhys());
367  m_fields[i]->SetPhysState(true);
368  transformed[i] = true;
369  }
370  }
371  Checkpoint_Output(nchk++);
372  for(i = 0; i < nfields; i++)
373  {
374  if (transformed[i])
375  {
376  m_fields[i]->SetWaveSpace(true);
377  m_fields[i]->HomogeneousFwdTrans(
378  m_fields[i]->GetPhys(),
379  m_fields[i]->UpdatePhys());
380  m_fields[i]->SetPhysState(false);
381  }
382  }
383  }
384  else
385  {
386  Checkpoint_Output(nchk++);
387  }
388  doCheckTime = false;
389  }
390 
391  // Step advance
392  ++step;
393  }
394 
395  // Print out summary statistics
396  if (m_session->GetComm()->GetRank() == 0)
397  {
398  if (m_cflSafetyFactor > 0.0)
399  {
400  cout << "CFL safety factor : " << m_cflSafetyFactor << endl
401  << "CFL time-step : " << m_timestep << endl;
402  }
403 
404  if (m_session->GetSolverInfo("Driver") != "SteadyState")
405  {
406  cout << "Time-integration : " << intTime << "s" << endl;
407  }
408  }
409 
410  // If homogeneous, transform back into physical space if necessary.
412  {
413  for(i = 0; i < nfields; i++)
414  {
415  if (m_fields[i]->GetWaveSpace())
416  {
417  m_fields[i]->SetWaveSpace(false);
418  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
419  m_fields[i]->UpdatePhys());
420  m_fields[i]->SetPhysState(true);
421  }
422  }
423  }
424  else
425  {
426  for(i = 0; i < nvariables; ++i)
427  {
428  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
429  m_fields[m_intVariables[i]]->SetPhysState(true);
430  }
431  }
432 
433  // Finalise filters
434  for (x = m_filters.begin(); x != m_filters.end(); ++x)
435  {
436  (*x)->Finalise(m_fields, m_time);
437  }
438 
439  // Print for 1D problems
440  if(m_spacedim == 1)
441  {
442  v_AppendOutput1D(fields);
443  }
444  }
445 
446  /**
447  * @brief Sets the initial conditions.
448  */
450  {
454  }
455 
456  /**
457  * @brief Prints a summary with some information regards the
458  * time-stepping.
459  */
461  {
463  AddSummaryItem(s, "Advection",
464  m_explicitAdvection ? "explicit" : "implicit");
465  AddSummaryItem(s, "Diffusion",
466  m_explicitDiffusion ? "explicit" : "implicit");
467 
468  if (m_session->GetSolverInfo("EQTYPE")
469  == "SteadyAdvectionDiffusionReaction")
470  {
471  AddSummaryItem(s, "Reaction",
472  m_explicitReaction ? "explicit" : "implicit");
473  }
474 
475  AddSummaryItem(s, "Time Step", m_timestep);
476  AddSummaryItem(s, "No. of Steps", m_steps);
477  AddSummaryItem(s, "Checkpoints (steps)", m_checksteps);
478  AddSummaryItem(s, "Integration Type",
480  m_intScheme->GetIntegrationMethod()]);
481  }
482 
483  /**
484  * Stores the solution in a file for 1D problems only. This method has
485  * been implemented to facilitate the post-processing for 1D problems.
486  */
488  Array<OneD, Array<OneD, NekDouble> > &solution1D)
489  {
490  // Coordinates of the quadrature points in the real physical space
494  m_fields[0]->GetCoords(x, y, z);
495 
496  // Print out the solution in a txt file
497  ofstream outfile;
498  outfile.open("solution1D.txt");
499  for(int i = 0; i < GetNpoints(); i++)
500  {
501  outfile << scientific << setw (17) << setprecision(16) << x[i]
502  << " " << solution1D[0][i] << endl;
503  }
504  outfile << endl << endl;
505  outfile.close();
506  }
507 
509  Array<OneD, Array<OneD, NekDouble> > &physfield,
510  Array<OneD, Array<OneD, NekDouble> > &numflux)
511  {
512  ASSERTL0(false,
513  "This function is not implemented for this equation.");
514  }
515 
517  Array<OneD, Array<OneD, NekDouble> > &physfield,
518  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
519  Array<OneD, Array<OneD, NekDouble> > &numfluxY )
520  {
521  ASSERTL0(false,
522  "This function is not implemented for this equation.");
523  }
524 
526  const Array<OneD, Array<OneD, NekDouble> > &ufield,
528  {
529  int i, j;
530  int nTraceNumPoints = GetTraceNpoints();
531  int nvariables = m_fields.num_elements();
532  int nqvar = uflux.num_elements();
533 
534  Array<OneD, NekDouble > Fwd (nTraceNumPoints);
535  Array<OneD, NekDouble > Bwd (nTraceNumPoints);
536  Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
537  Array<OneD, NekDouble > fluxtemp(nTraceNumPoints, 0.0);
538 
539  // Get the sign of (v \cdot n), v = an arbitrary vector
540 
541  // Evaulate upwind flux:
542  // uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
543  for (j = 0; j < nqvar; ++j)
544  {
545  for (i = 0; i < nvariables ; ++i)
546  {
547  // Compute Fwd and Bwd value of ufield of i direction
548  m_fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
549 
550  // if Vn >= 0, flux = uFwd, i.e.,
551  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uFwd
552  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uFwd
553 
554  // else if Vn < 0, flux = uBwd, i.e.,
555  // edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uBwd
556  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uBwd
557 
558  m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
559  Fwd, Bwd, fluxtemp);
560 
561  // Imposing weak boundary condition with flux
562  // if Vn >= 0, uflux = uBwd at Neumann, i.e.,
563  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uBwd
564  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uBwd
565 
566  // if Vn >= 0, uflux = uFwd at Neumann, i.e.,
567  // edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uFwd
568  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uFwd
569 
570  if(m_fields[0]->GetBndCondExpansions().num_elements())
571  {
572  WeakPenaltyforScalar(i, ufield[i], fluxtemp);
573  }
574 
575  // if Vn >= 0, flux = uFwd*(tan_{\xi}^- \cdot \vec{n}),
576  // i.e,
577  // edge::eForward, uFwd \(\tan_{\xi}^Fwd \cdot \vec{n})
578  // edge::eBackward, uFwd \(\tan_{\xi}^Bwd \cdot \vec{n})
579 
580  // else if Vn < 0, flux = uBwd*(tan_{\xi}^- \cdot \vec{n}),
581  // i.e,
582  // edge::eForward, uBwd \(\tan_{\xi}^Fwd \cdot \vec{n})
583  // edge::eBackward, uBwd \(\tan_{\xi}^Bwd \cdot \vec{n})
584 
585  Vmath::Vmul(nTraceNumPoints,
586  m_traceNormals[j], 1,
587  fluxtemp, 1,
588  uflux[j][i], 1);
589  }
590  }
591  }
592 
593 
594 
596  const Array<OneD, Array<OneD, NekDouble> > &ufield,
599  {
600  int nTraceNumPoints = GetTraceNpoints();
601  int nvariables = m_fields.num_elements();
602  int nqvar = qfield.num_elements();
603 
604  NekDouble C11 = 1.0;
605  Array<OneD, NekDouble > Fwd(nTraceNumPoints);
606  Array<OneD, NekDouble > Bwd(nTraceNumPoints);
607  Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
608 
609  Array<OneD, NekDouble > qFwd (nTraceNumPoints);
610  Array<OneD, NekDouble > qBwd (nTraceNumPoints);
611  Array<OneD, NekDouble > qfluxtemp(nTraceNumPoints, 0.0);
612 
613  Array<OneD, NekDouble > uterm(nTraceNumPoints);
614 
615  // Evaulate upwind flux:
616  // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
617  for (int i = 0; i < nvariables; ++i)
618  {
619  qflux[i] = Array<OneD, NekDouble> (nTraceNumPoints, 0.0);
620  for (int j = 0; j < nqvar; ++j)
621  {
622  // Compute Fwd and Bwd value of ufield of jth direction
623  m_fields[i]->GetFwdBwdTracePhys(qfield[j][i],qFwd,qBwd);
624 
625  // if Vn >= 0, flux = uFwd, i.e.,
626  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick
627  // qflux = qBwd = q+
628  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick
629  // qflux = qBwd = q-
630 
631  // else if Vn < 0, flux = uBwd, i.e.,
632  // edge::eForward, if V*n<0 <=> V*n_F<0, pick
633  // qflux = qFwd = q-
634  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick
635  // qflux = qFwd = q+
636 
637  m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
638  qBwd, qFwd,
639  qfluxtemp);
640 
641  Vmath::Vmul(nTraceNumPoints,
642  m_traceNormals[j], 1,
643  qfluxtemp, 1,
644  qfluxtemp, 1);
645 
646  // Generate Stability term = - C11 ( u- - u+ )
647  m_fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
648 
649  Vmath::Vsub(nTraceNumPoints,
650  Fwd, 1, Bwd, 1,
651  uterm, 1);
652 
653  Vmath::Smul(nTraceNumPoints,
654  -1.0 * C11, uterm, 1,
655  uterm, 1);
656 
657  // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
658  Vmath::Vadd(nTraceNumPoints,
659  uterm, 1,
660  qfluxtemp, 1,
661  qfluxtemp, 1);
662 
663  // Imposing weak boundary condition with flux
664  if (m_fields[0]->GetBndCondExpansions().num_elements())
665  {
666  WeakPenaltyforVector(i, j,
667  qfield[j][i],
668  qfluxtemp,
669  C11);
670  }
671 
672  // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
673  // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
674  // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
675  Vmath::Vadd(nTraceNumPoints,
676  qfluxtemp, 1,
677  qflux[i], 1,
678  qflux[i], 1);
679  }
680  }
681  }
682 
684  {
685  if (m_session->DefinesFunction("InitialConditions"))
686  {
687  for (int i = 0; i < m_fields.num_elements(); ++i)
688  {
690 
691  vType = m_session->GetFunctionType(
692  "InitialConditions", m_session->GetVariable(i));
693 
694  if (vType == LibUtilities::eFunctionTypeFile)
695  {
696  std::string filename
697  = m_session->GetFunctionFilename(
698  "InitialConditions", m_session->GetVariable(i));
699 
700  fs::path pfilename(filename);
701 
702  // redefine path for parallel file which is in directory
703  if(fs::is_directory(pfilename))
704  {
705  fs::path metafile("Info.xml");
706  fs::path fullpath = pfilename / metafile;
707  filename = LibUtilities::PortablePath(fullpath);
708  }
709  m_fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
710 
711  // check to see if time defined
712  if (m_fieldMetaDataMap !=
714  {
716 
717  iter = m_fieldMetaDataMap.find("Time");
718  if (iter != m_fieldMetaDataMap.end())
719  {
720  time = boost::lexical_cast<NekDouble>(
721  iter->second);
722  }
723  }
724 
725  break;
726  }
727  }
728  }
729  }
730 
732  const int var,
733  const Array<OneD, const NekDouble> &physfield,
734  Array<OneD, NekDouble> &penaltyflux,
735  NekDouble time)
736  {
737  int i, e, npoints, id1, id2;
738 
739  // Number of boundary regions
740  int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
741  int Nfps, numBDEdge;
742  int nTraceNumPoints = GetTraceNpoints();
743  int cnt = 0;
744 
745  Array<OneD, NekDouble > uplus(nTraceNumPoints);
746 
747  m_fields[var]->ExtractTracePhys(physfield, uplus);
748  for (i = 0; i < nbnd; ++i)
749  {
750  // Number of boundary expansion related to that region
751  numBDEdge = m_fields[var]->
752  GetBndCondExpansions()[i]->GetExpSize();
753 
754  // Evaluate boundary values g_D or g_N from input files
756  m_session->GetFunction("InitialConditions", 0);
757 
758  npoints = m_fields[var]->
759  GetBndCondExpansions()[i]->GetNpoints();
760 
761  Array<OneD,NekDouble> BDphysics(npoints);
762  Array<OneD,NekDouble> x0(npoints,0.0);
763  Array<OneD,NekDouble> x1(npoints,0.0);
764  Array<OneD,NekDouble> x2(npoints,0.0);
765 
766  m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
767  ifunc->Evaluate(x0,x1,x2,time,BDphysics);
768 
769  // Weakly impose boundary conditions by modifying flux values
770  for (e = 0; e < numBDEdge ; ++e)
771  {
772  // Number of points on the expansion
773  Nfps = m_fields[var]->
774  GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
775 
776  id1 = m_fields[var]->
777  GetBndCondExpansions()[i]->GetPhys_Offset(e);
778 
779  id2 = m_fields[0]->GetTrace()->
780  GetPhys_Offset(m_fields[0]->GetTraceMap()->
781  GetBndCondTraceToGlobalTraceMap(cnt++));
782 
783  // For Dirichlet boundary condition: uflux = g_D
784  if (m_fields[var]->GetBndConditions()[i]->
785  GetBoundaryConditionType() == SpatialDomains::eDirichlet)
786  {
787  Vmath::Vcopy(Nfps,
788  &BDphysics[id1], 1,
789  &penaltyflux[id2], 1);
790  }
791  // For Neumann boundary condition: uflux = u+
792  else if ((m_fields[var]->GetBndConditions()[i])->
793  GetBoundaryConditionType() == SpatialDomains::eNeumann)
794  {
795  Vmath::Vcopy(Nfps,
796  &uplus[id2], 1,
797  &penaltyflux[id2], 1);
798  }
799  }
800  }
801  }
802 
803  /**
804  * Diffusion: Imposing weak boundary condition for q with flux
805  * uflux = g_D on Dirichlet boundary condition
806  * uflux = u_Fwd on Neumann boundary condition
807  */
809  const int var,
810  const int dir,
811  const Array<OneD, const NekDouble> &physfield,
812  Array<OneD, NekDouble> &penaltyflux,
813  NekDouble C11,
814  NekDouble time)
815  {
816  int i, e, npoints, id1, id2;
817  int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
818  int numBDEdge, Nfps;
819  int nTraceNumPoints = GetTraceNpoints();
820  Array<OneD, NekDouble > uterm(nTraceNumPoints);
821  Array<OneD, NekDouble > qtemp(nTraceNumPoints);
822  int cnt = 0;
823 
824  m_fields[var]->ExtractTracePhys(physfield,qtemp);
825 
826  for (i = 0; i < nbnd; ++i)
827  {
828  numBDEdge = m_fields[var]->
829  GetBndCondExpansions()[i]->GetExpSize();
830 
831  // Evaluate boundary values g_D or g_N from input files
833  m_session->GetFunction("InitialConditions", 0);
834 
835  npoints = m_fields[var]->
836  GetBndCondExpansions()[i]->GetNpoints();
837 
838  Array<OneD,NekDouble> BDphysics(npoints);
839  Array<OneD,NekDouble> x0(npoints,0.0);
840  Array<OneD,NekDouble> x1(npoints,0.0);
841  Array<OneD,NekDouble> x2(npoints,0.0);
842 
843  m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
844  ifunc->Evaluate(x0,x1,x2,time,BDphysics);
845 
846  // Weakly impose boundary conditions by modifying flux values
847  for (e = 0; e < numBDEdge ; ++e)
848  {
849  Nfps = m_fields[var]->
850  GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
851 
852  id1 = m_fields[var]->
853  GetBndCondExpansions()[i]->GetPhys_Offset(e);
854 
855  id2 = m_fields[0]->GetTrace()->
856  GetPhys_Offset(m_fields[0]->GetTraceMap()->
857  GetBndCondTraceToGlobalTraceMap(cnt++));
858 
859  // For Dirichlet boundary condition:
860  //qflux = q+ - C_11 (u+ - g_D) (nx, ny)
861  if(m_fields[var]->GetBndConditions()[i]->
862  GetBoundaryConditionType() == SpatialDomains::eDirichlet)
863  {
864  Vmath::Vmul(Nfps,
865  &m_traceNormals[dir][id2], 1,
866  &qtemp[id2], 1,
867  &penaltyflux[id2], 1);
868  }
869  // For Neumann boundary condition: qflux = g_N
870  else if((m_fields[var]->GetBndConditions()[i])->
871  GetBoundaryConditionType() == SpatialDomains::eNeumann)
872  {
873  Vmath::Vmul(Nfps,
874  &m_traceNormals[dir][id2], 1,
875  &BDphysics[id1], 1,
876  &penaltyflux[id2], 1);
877  }
878  }
879  }
880  }
881 
882  /**
883  * @brief Return the timestep to be used for the next step in the
884  * time-marching loop.
885  *
886  * This function can be overloaded to facilitate solver which utilise a
887  * CFL (or other) parameter to determine a maximum timestep under which
888  * the problem remains stable.
889  */
891  const Array<OneD, const Array<OneD, NekDouble> > &inarray)
892  {
893  return v_GetTimeStep(inarray);
894  }
895 
896  /**
897  * @brief Return the timestep to be used for the next step in the
898  * time-marching loop.
899  *
900  * @see UnsteadySystem::GetTimeStep
901  */
903  const Array<OneD, const Array<OneD, NekDouble> > &inarray)
904  {
905  ASSERTL0(false, "Not defined for this class");
906  return 0.0;
907  }
908 
910  {
911  return false;
912  }
913 
915  {
916  return false;
917  }
918 
920  {
921  return false;
922  }
923  }
924 }
virtual SOLVER_UTILS_EXPORT bool v_SteadyStateCheck(int step)
virtual SOLVER_UTILS_EXPORT void v_DoSolve()
Solves an unsteady problem.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &l)
Virtual function for generating summary information.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
bool m_explicitReaction
Indicates if explicit or implicit treatment of reaction is used.
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector(const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
SOLVER_UTILS_EXPORT void CheckForRestartTime(NekDouble &time)
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate(int step)
A base class for describing how to solve specific equations.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
Adams-Bashforth Forward multi-step scheme of order 2.
NekDouble m_time
Current time of simulation.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
Runge-Kutta multi-stage scheme 4th order explicit (old name)
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
NekDouble m_timestep
Time step size.
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:47
void WeakPenaltyforVector(const int var, const int dir, const Array< OneD, const NekDouble > &physfield, Array< OneD, NekDouble > &penaltyflux, NekDouble C11, NekDouble time=0.0)
Classical RungeKutta4 method (new name for eClassicalRungeKutta4)
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D(Array< OneD, Array< OneD, NekDouble > > &solution1D)
Print the solution at each solution point in a txt file.
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep(const Array< OneD, const Array< OneD, NekDouble > > &inarray)
Return the timestep to be used for the next step in the time-marching loop.
NekDouble m_checktime
Time between checkpoints.
const char *const TimeIntegrationMethodMap[]
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator()
Get the maximum timestep estimator for cfl control.
Nonlinear SSP RungeKutta3 explicit.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:50
int m_checksteps
Number of steps between checkpoints.
NekDouble m_fintime
Finish time of the simulation.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Initialisation object for EquationSystem.
int m_steps
Number of steps to take.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
Classical RungeKutta2 method (new name for eMidpoint)
void Stop()
Definition: Timer.cpp:62
static const NekDouble kNekZeroTol
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
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem()
Destructor.
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate(int step)
Adams-Bashforth Forward multi-step scheme of order 1.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:50
int Nnan(int n, const T *x, const int incx)
Return number of NaN elements of x.
Definition: Vmath.cpp:869
std::vector< std::pair< std::string, FilterParams > > FilterMap
Definition: SessionReader.h:66
Nonlinear SSP RungeKutta2 explicit (surrogate for eRungeKutta2_ImprovedEuler)
Improved RungeKutta2 explicit (old name meaning Heun's method)
int m_spacedim
Spatial dimension (>= expansion dim).
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
Definition: FileSystem.cpp:41
double NekDouble
TimeIntegrationWrapperFactory & GetTimeIntegrationWrapperFactory()
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep()
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar(const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
boost::shared_ptr< Equation > EquationSharedPtr
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises UnsteadySystem class members.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
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
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
SOLVER_UTILS_EXPORT int GetPhys_Offset(int n)
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
LibUtilities::FieldIOSharedPtr m_fld
Field input/output.
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:42
std::vector< FilterSharedPtr > m_filters
int m_infosteps
Number of time steps between outputting status information.
SOLVER_UTILS_EXPORT void SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Initialise the data in the dependent fields.
void Start()
Definition: Timer.cpp:51
void WeakPenaltyforScalar(const int var, const Array< OneD, const NekDouble > &physfield, Array< OneD, NekDouble > &penaltyflux, NekDouble time=0.0)
midpoint method (old name)
NekDouble TimePerTest(unsigned int n)
Returns amount of seconds per iteration in a test with n iterations.
Definition: Timer.cpp:108
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:66
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
virtual SOLVER_UTILS_EXPORT void v_DoInitialise()
Sets up initial conditions.
enum HomogeneousType m_HomogeneousType
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