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  {
136  {
137  TimeStability = 2.784;
138  break;
139  }
143  {
144  TimeStability = 2.0;
145  break;
146  }
148  {
149  TimeStability = 1.0;
150  break;
151  }
152  default:
153  {
154  ASSERTL0(
155  false,
156  "No CFL control implementation for this time"
157  "integration scheme");
158  }
159  }
160  return TimeStability;
161  }
162 
163  /**
164  * @brief Initialises the time integration scheme (as specified in the
165  * session file), and perform the time integration.
166  */
168  {
169  ASSERTL0(m_intScheme != 0, "No time integration scheme.");
170 
171  int i, nchk = 1;
172  int nvariables = 0;
173  int nfields = m_fields.num_elements();
174 
175  if (m_intVariables.empty())
176  {
177  for (i = 0; i < nfields; ++i)
178  {
179  m_intVariables.push_back(i);
180  }
181  nvariables = nfields;
182  }
183  else
184  {
185  nvariables = m_intVariables.size();
186  }
187 
188  // Integrate in wave-space if using homogeneous1D
190  {
191  for(i = 0; i < nfields; ++i)
192  {
193  m_fields[i]->HomogeneousFwdTrans(m_fields[i]->GetPhys(),
194  m_fields[i]->UpdatePhys());
195  m_fields[i]->SetWaveSpace(true);
196  m_fields[i]->SetPhysState(false);
197  }
198  }
199 
200  // Set up wrapper to fields data storage.
201  Array<OneD, Array<OneD, NekDouble> > fields(nvariables);
202  Array<OneD, Array<OneD, NekDouble> > tmp (nvariables);
203 
204  // Order storage to list time-integrated fields first.
205  for(i = 0; i < nvariables; ++i)
206  {
207  fields[i] = m_fields[m_intVariables[i]]->GetPhys();
208  m_fields[m_intVariables[i]]->SetPhysState(false);
209  }
210 
211  // Initialise time integration scheme
212  m_intSoln = m_intScheme->InitializeScheme(
213  m_timestep, fields, m_time, m_ode);
214 
215  // Initialise filters
217  for (x = m_filters.begin(); x != m_filters.end(); ++x)
218  {
219  (*x)->Initialise(m_fields, m_time);
220  }
221 
222  // Ensure that there is no conflict of parameters
223  if(m_cflSafetyFactor > 0.0)
224  {
225  // Check final condition
226  ASSERTL0(m_fintime == 0.0 || m_steps == 0,
227  "Final condition not unique: "
228  "fintime > 0.0 and Nsteps > 0");
229 
230  // Check timestep condition
231  ASSERTL0(m_timestep == 0.0,
232  "Timestep not unique: timestep > 0.0 & CFL > 0.0");
233  }
234 
235  // Check uniqueness of checkpoint output
236  ASSERTL0((m_checktime == 0.0 && m_checksteps == 0) ||
237  (m_checktime > 0.0 && m_checksteps == 0) ||
238  (m_checktime == 0.0 && m_checksteps > 0),
239  "Only one of IO_CheckTime and IO_CheckSteps "
240  "should be set!");
241 
242  Timer timer;
243  bool doCheckTime = false;
244  int step = 0;
245  NekDouble intTime = 0.0;
246  NekDouble lastCheckTime = 0.0;
247  NekDouble cpuTime = 0.0;
248  NekDouble elapsed = 0.0;
249 
250  while (step < m_steps ||
252  {
253  if (m_cflSafetyFactor)
254  {
255  m_timestep = GetTimeStep(fields);
256 
257  // Ensure that the final timestep finishes at the final
258  // time, or at a prescribed IO_CheckTime.
259  if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
260  {
262  }
263  else if (m_checktime &&
264  m_time + m_timestep - lastCheckTime >= m_checktime)
265  {
266  lastCheckTime += m_checktime;
267  m_timestep = lastCheckTime - m_time;
268  doCheckTime = true;
269  }
270  }
271 
272  // Perform any solver-specific pre-integration steps
273  if (v_PreIntegrate(step))
274  {
275  break;
276  }
277 
278  timer.Start();
279  fields = m_intScheme->TimeIntegrate(
280  step, m_timestep, m_intSoln, m_ode);
281  timer.Stop();
282 
283  m_time += m_timestep;
284  elapsed = timer.TimePerTest(1);
285  intTime += elapsed;
286  cpuTime += elapsed;
287 
288  // Write out status information
289  if (m_session->GetComm()->GetRank() == 0 &&
290  !((step+1) % m_infosteps))
291  {
292  cout << "Steps: " << setw(8) << left << step+1 << " "
293  << "Time: " << setw(12) << left << m_time;
294 
295  if (m_cflSafetyFactor)
296  {
297  cout << " Time-step: " << setw(12)
298  << left << m_timestep;
299  }
300 
301  stringstream ss;
302  ss << cpuTime << "s";
303  cout << " CPU Time: " << setw(8) << left
304  << ss.str() << endl;
305  cpuTime = 0.0;
306  }
307 
308  // Transform data into coefficient space
309  for (i = 0; i < nvariables; ++i)
310  {
311  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
312  m_fields[m_intVariables[i]]->FwdTrans_IterPerExp(
313  fields[i],
314  m_fields[m_intVariables[i]]->UpdateCoeffs());
315  m_fields[m_intVariables[i]]->SetPhysState(false);
316  }
317 
318  // Perform any solver-specific post-integration steps
319  if (v_PostIntegrate(step))
320  {
321  break;
322  }
323 
324  // Update filters
326  for (x = m_filters.begin(); x != m_filters.end(); ++x)
327  {
328  (*x)->Update(m_fields, m_time);
329  }
330 
331  // Write out checkpoint files
332  if ((m_checksteps && step && !((step + 1) % m_checksteps)) ||
333  doCheckTime)
334  {
336  {
337  vector<bool> transformed(nfields, false);
338  for(i = 0; i < nfields; i++)
339  {
340  if (m_fields[i]->GetWaveSpace())
341  {
342  m_fields[i]->SetWaveSpace(false);
343  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
344  m_fields[i]->UpdatePhys());
345  m_fields[i]->SetPhysState(true);
346  transformed[i] = true;
347  }
348  }
349  Checkpoint_Output(nchk++);
350  for(i = 0; i < nfields; i++)
351  {
352  if (transformed[i])
353  {
354  m_fields[i]->SetWaveSpace(true);
355  m_fields[i]->HomogeneousFwdTrans(
356  m_fields[i]->GetPhys(),
357  m_fields[i]->UpdatePhys());
358  m_fields[i]->SetPhysState(false);
359  }
360  }
361  }
362  else
363  {
364  Checkpoint_Output(nchk++);
365  }
366  doCheckTime = false;
367  }
368 
369  // Step advance
370  ++step;
371  }
372 
373  // Print out summary statistics
374  if (m_session->GetComm()->GetRank() == 0)
375  {
376  if (m_cflSafetyFactor > 0.0)
377  {
378  cout << "CFL safety factor : " << m_cflSafetyFactor << endl
379  << "CFL time-step : " << m_timestep << endl;
380  }
381 
382  if (m_session->GetSolverInfo("Driver") != "SteadyState")
383  {
384  cout << "Time-integration : " << intTime << "s" << endl;
385  }
386  }
387 
388  // If homogeneous, transform back into physical space if necessary.
390  {
391  for(i = 0; i < nfields; i++)
392  {
393  if (m_fields[i]->GetWaveSpace())
394  {
395  m_fields[i]->SetWaveSpace(false);
396  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
397  m_fields[i]->UpdatePhys());
398  m_fields[i]->SetPhysState(true);
399  }
400  }
401  }
402  else
403  {
404  for(i = 0; i < nvariables; ++i)
405  {
406  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
407  m_fields[m_intVariables[i]]->SetPhysState(true);
408  }
409  }
410 
411  // Finalise filters
412  for (x = m_filters.begin(); x != m_filters.end(); ++x)
413  {
414  (*x)->Finalise(m_fields, m_time);
415  }
416 
417  // Print for 1D problems
418  if(m_spacedim == 1)
419  {
420  v_AppendOutput1D(fields);
421  }
422  }
423 
424  /**
425  * @brief Sets the initial conditions.
426  */
428  {
432  }
433 
434  /**
435  * @brief Prints a summary with some information regards the
436  * time-stepping.
437  */
439  {
441  AddSummaryItem(s, "Advection",
442  m_explicitAdvection ? "explicit" : "implicit");
443  AddSummaryItem(s, "Diffusion",
444  m_explicitDiffusion ? "explicit" : "implicit");
445 
446  if (m_session->GetSolverInfo("EQTYPE")
447  == "SteadyAdvectionDiffusionReaction")
448  {
449  AddSummaryItem(s, "Reaction",
450  m_explicitReaction ? "explicit" : "implicit");
451  }
452 
453  AddSummaryItem(s, "Time Step", m_timestep);
454  AddSummaryItem(s, "No. of Steps", m_steps);
455  AddSummaryItem(s, "Checkpoints (steps)", m_checksteps);
456  AddSummaryItem(s, "Integration Type",
458  m_intScheme->GetIntegrationMethod()]);
459  }
460 
461  /**
462  * Stores the solution in a file for 1D problems only. This method has
463  * been implemented to facilitate the post-processing for 1D problems.
464  */
466  Array<OneD, Array<OneD, NekDouble> > &solution1D)
467  {
468  // Coordinates of the quadrature points in the real physical space
472  m_fields[0]->GetCoords(x, y, z);
473 
474  // Print out the solution in a txt file
475  ofstream outfile;
476  outfile.open("solution1D.txt");
477  for(int i = 0; i < GetNpoints(); i++)
478  {
479  outfile << scientific << setw (17) << setprecision(16) << x[i]
480  << " " << solution1D[0][i] << endl;
481  }
482  outfile << endl << endl;
483  outfile.close();
484  }
485 
487  Array<OneD, Array<OneD, NekDouble> > &physfield,
488  Array<OneD, Array<OneD, NekDouble> > &numflux)
489  {
490  ASSERTL0(false,
491  "This function is not implemented for this equation.");
492  }
493 
495  Array<OneD, Array<OneD, NekDouble> > &physfield,
496  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
497  Array<OneD, Array<OneD, NekDouble> > &numfluxY )
498  {
499  ASSERTL0(false,
500  "This function is not implemented for this equation.");
501  }
502 
504  const Array<OneD, Array<OneD, NekDouble> > &ufield,
506  {
507  int i, j;
508  int nTraceNumPoints = GetTraceNpoints();
509  int nvariables = m_fields.num_elements();
510  int nqvar = uflux.num_elements();
511 
512  Array<OneD, NekDouble > Fwd (nTraceNumPoints);
513  Array<OneD, NekDouble > Bwd (nTraceNumPoints);
514  Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
515  Array<OneD, NekDouble > fluxtemp(nTraceNumPoints, 0.0);
516 
517  // Get the sign of (v \cdot n), v = an arbitrary vector
518 
519  // Evaulate upwind flux:
520  // uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
521  for (j = 0; j < nqvar; ++j)
522  {
523  for (i = 0; i < nvariables ; ++i)
524  {
525  // Compute Fwd and Bwd value of ufield of i direction
526  m_fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
527 
528  // if Vn >= 0, flux = uFwd, i.e.,
529  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uFwd
530  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uFwd
531 
532  // else if Vn < 0, flux = uBwd, i.e.,
533  // edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uBwd
534  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uBwd
535 
536  m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
537  Fwd, Bwd, fluxtemp);
538 
539  // Imposing weak boundary condition with flux
540  // if Vn >= 0, uflux = uBwd at Neumann, i.e.,
541  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uBwd
542  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uBwd
543 
544  // if Vn >= 0, uflux = uFwd at Neumann, i.e.,
545  // edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uFwd
546  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uFwd
547 
548  if(m_fields[0]->GetBndCondExpansions().num_elements())
549  {
550  WeakPenaltyforScalar(i, ufield[i], fluxtemp);
551  }
552 
553  // if Vn >= 0, flux = uFwd*(tan_{\xi}^- \cdot \vec{n}),
554  // i.e,
555  // edge::eForward, uFwd \(\tan_{\xi}^Fwd \cdot \vec{n})
556  // edge::eBackward, uFwd \(\tan_{\xi}^Bwd \cdot \vec{n})
557 
558  // else if Vn < 0, flux = uBwd*(tan_{\xi}^- \cdot \vec{n}),
559  // i.e,
560  // edge::eForward, uBwd \(\tan_{\xi}^Fwd \cdot \vec{n})
561  // edge::eBackward, uBwd \(\tan_{\xi}^Bwd \cdot \vec{n})
562 
563  Vmath::Vmul(nTraceNumPoints,
564  m_traceNormals[j], 1,
565  fluxtemp, 1,
566  uflux[j][i], 1);
567  }
568  }
569  }
570 
571 
572 
574  const Array<OneD, Array<OneD, NekDouble> > &ufield,
577  {
578  int nTraceNumPoints = GetTraceNpoints();
579  int nvariables = m_fields.num_elements();
580  int nqvar = qfield.num_elements();
581 
582  NekDouble C11 = 1.0;
583  Array<OneD, NekDouble > Fwd(nTraceNumPoints);
584  Array<OneD, NekDouble > Bwd(nTraceNumPoints);
585  Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
586 
587  Array<OneD, NekDouble > qFwd (nTraceNumPoints);
588  Array<OneD, NekDouble > qBwd (nTraceNumPoints);
589  Array<OneD, NekDouble > qfluxtemp(nTraceNumPoints, 0.0);
590 
591  Array<OneD, NekDouble > uterm(nTraceNumPoints);
592 
593  // Evaulate upwind flux:
594  // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
595  for (int i = 0; i < nvariables; ++i)
596  {
597  qflux[i] = Array<OneD, NekDouble> (nTraceNumPoints, 0.0);
598  for (int j = 0; j < nqvar; ++j)
599  {
600  // Compute Fwd and Bwd value of ufield of jth direction
601  m_fields[i]->GetFwdBwdTracePhys(qfield[j][i],qFwd,qBwd);
602 
603  // if Vn >= 0, flux = uFwd, i.e.,
604  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick
605  // qflux = qBwd = q+
606  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick
607  // qflux = qBwd = q-
608 
609  // else if Vn < 0, flux = uBwd, i.e.,
610  // edge::eForward, if V*n<0 <=> V*n_F<0, pick
611  // qflux = qFwd = q-
612  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick
613  // qflux = qFwd = q+
614 
615  m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
616  qBwd, qFwd,
617  qfluxtemp);
618 
619  Vmath::Vmul(nTraceNumPoints,
620  m_traceNormals[j], 1,
621  qfluxtemp, 1,
622  qfluxtemp, 1);
623 
624  // Generate Stability term = - C11 ( u- - u+ )
625  m_fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
626 
627  Vmath::Vsub(nTraceNumPoints,
628  Fwd, 1, Bwd, 1,
629  uterm, 1);
630 
631  Vmath::Smul(nTraceNumPoints,
632  -1.0 * C11, uterm, 1,
633  uterm, 1);
634 
635  // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
636  Vmath::Vadd(nTraceNumPoints,
637  uterm, 1,
638  qfluxtemp, 1,
639  qfluxtemp, 1);
640 
641  // Imposing weak boundary condition with flux
642  if (m_fields[0]->GetBndCondExpansions().num_elements())
643  {
644  WeakPenaltyforVector(i, j,
645  qfield[j][i],
646  qfluxtemp,
647  C11);
648  }
649 
650  // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
651  // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
652  // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
653  Vmath::Vadd(nTraceNumPoints,
654  qfluxtemp, 1,
655  qflux[i], 1,
656  qflux[i], 1);
657  }
658  }
659  }
660 
662  {
663  if (m_session->DefinesFunction("InitialConditions"))
664  {
665  for (int i = 0; i < m_fields.num_elements(); ++i)
666  {
668 
669  vType = m_session->GetFunctionType(
670  "InitialConditions", m_session->GetVariable(i));
671 
672  if (vType == LibUtilities::eFunctionTypeFile)
673  {
674  std::string filename
675  = m_session->GetFunctionFilename(
676  "InitialConditions", m_session->GetVariable(i));
677 
678  fs::path pfilename(filename);
679 
680  // redefine path for parallel file which is in directory
681  if(fs::is_directory(pfilename))
682  {
683  fs::path metafile("Info.xml");
684  fs::path fullpath = pfilename / metafile;
685  filename = LibUtilities::PortablePath(fullpath);
686  }
687  m_fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
688 
689  // check to see if time defined
690  if (m_fieldMetaDataMap !=
692  {
694 
695  iter = m_fieldMetaDataMap.find("Time");
696  if (iter != m_fieldMetaDataMap.end())
697  {
698  time = boost::lexical_cast<NekDouble>(
699  iter->second);
700  }
701  }
702 
703  break;
704  }
705  }
706  }
707  }
708 
710  const int var,
711  const Array<OneD, const NekDouble> &physfield,
712  Array<OneD, NekDouble> &penaltyflux,
713  NekDouble time)
714  {
715  int i, e, npoints, id1, id2;
716 
717  // Number of boundary regions
718  int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
719  int Nfps, numBDEdge;
720  int nTraceNumPoints = GetTraceNpoints();
721  int cnt = 0;
722 
723  Array<OneD, NekDouble > uplus(nTraceNumPoints);
724 
725  m_fields[var]->ExtractTracePhys(physfield, uplus);
726  for (i = 0; i < nbnd; ++i)
727  {
728  // Number of boundary expansion related to that region
729  numBDEdge = m_fields[var]->
730  GetBndCondExpansions()[i]->GetExpSize();
731 
732  // Evaluate boundary values g_D or g_N from input files
734  m_session->GetFunction("InitialConditions", 0);
735 
736  npoints = m_fields[var]->
737  GetBndCondExpansions()[i]->GetNpoints();
738 
739  Array<OneD,NekDouble> BDphysics(npoints);
740  Array<OneD,NekDouble> x0(npoints,0.0);
741  Array<OneD,NekDouble> x1(npoints,0.0);
742  Array<OneD,NekDouble> x2(npoints,0.0);
743 
744  m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
745  ifunc->Evaluate(x0,x1,x2,time,BDphysics);
746 
747  // Weakly impose boundary conditions by modifying flux values
748  for (e = 0; e < numBDEdge ; ++e)
749  {
750  // Number of points on the expansion
751  Nfps = m_fields[var]->
752  GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
753 
754  id1 = m_fields[var]->
755  GetBndCondExpansions()[i]->GetPhys_Offset(e);
756 
757  id2 = m_fields[0]->GetTrace()->
758  GetPhys_Offset(m_fields[0]->GetTraceMap()->
759  GetBndCondTraceToGlobalTraceMap(cnt++));
760 
761  // For Dirichlet boundary condition: uflux = g_D
762  if (m_fields[var]->GetBndConditions()[i]->
763  GetBoundaryConditionType() == SpatialDomains::eDirichlet)
764  {
765  Vmath::Vcopy(Nfps,
766  &BDphysics[id1], 1,
767  &penaltyflux[id2], 1);
768  }
769  // For Neumann boundary condition: uflux = u+
770  else if ((m_fields[var]->GetBndConditions()[i])->
771  GetBoundaryConditionType() == SpatialDomains::eNeumann)
772  {
773  Vmath::Vcopy(Nfps,
774  &uplus[id2], 1,
775  &penaltyflux[id2], 1);
776  }
777  }
778  }
779  }
780 
781  /**
782  * Diffusion: Imposing weak boundary condition for q with flux
783  * uflux = g_D on Dirichlet boundary condition
784  * uflux = u_Fwd on Neumann boundary condition
785  */
787  const int var,
788  const int dir,
789  const Array<OneD, const NekDouble> &physfield,
790  Array<OneD, NekDouble> &penaltyflux,
791  NekDouble C11,
792  NekDouble time)
793  {
794  int i, e, npoints, id1, id2;
795  int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
796  int numBDEdge, Nfps;
797  int nTraceNumPoints = GetTraceNpoints();
798  Array<OneD, NekDouble > uterm(nTraceNumPoints);
799  Array<OneD, NekDouble > qtemp(nTraceNumPoints);
800  int cnt = 0;
801 
802  m_fields[var]->ExtractTracePhys(physfield,qtemp);
803 
804  for (i = 0; i < nbnd; ++i)
805  {
806  numBDEdge = m_fields[var]->
807  GetBndCondExpansions()[i]->GetExpSize();
808 
809  // Evaluate boundary values g_D or g_N from input files
811  m_session->GetFunction("InitialConditions", 0);
812 
813  npoints = m_fields[var]->
814  GetBndCondExpansions()[i]->GetNpoints();
815 
816  Array<OneD,NekDouble> BDphysics(npoints);
817  Array<OneD,NekDouble> x0(npoints,0.0);
818  Array<OneD,NekDouble> x1(npoints,0.0);
819  Array<OneD,NekDouble> x2(npoints,0.0);
820 
821  m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
822  ifunc->Evaluate(x0,x1,x2,time,BDphysics);
823 
824  // Weakly impose boundary conditions by modifying flux values
825  for (e = 0; e < numBDEdge ; ++e)
826  {
827  Nfps = m_fields[var]->
828  GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
829 
830  id1 = m_fields[var]->
831  GetBndCondExpansions()[i]->GetPhys_Offset(e);
832 
833  id2 = m_fields[0]->GetTrace()->
834  GetPhys_Offset(m_fields[0]->GetTraceMap()->
835  GetBndCondTraceToGlobalTraceMap(cnt++));
836 
837  // For Dirichlet boundary condition:
838  //qflux = q+ - C_11 (u+ - g_D) (nx, ny)
839  if(m_fields[var]->GetBndConditions()[i]->
840  GetBoundaryConditionType() == SpatialDomains::eDirichlet)
841  {
842  Vmath::Vmul(Nfps,
843  &m_traceNormals[dir][id2], 1,
844  &qtemp[id2], 1,
845  &penaltyflux[id2], 1);
846  }
847  // For Neumann boundary condition: qflux = g_N
848  else if((m_fields[var]->GetBndConditions()[i])->
849  GetBoundaryConditionType() == SpatialDomains::eNeumann)
850  {
851  Vmath::Vmul(Nfps,
852  &m_traceNormals[dir][id2], 1,
853  &BDphysics[id1], 1,
854  &penaltyflux[id2], 1);
855  }
856  }
857  }
858  }
859 
860  /**
861  * @brief Return the timestep to be used for the next step in the
862  * time-marching loop.
863  *
864  * This function can be overloaded to facilitate solver which utilise a
865  * CFL (or other) parameter to determine a maximum timestep under which
866  * the problem remains stable.
867  */
869  const Array<OneD, const Array<OneD, NekDouble> > &inarray)
870  {
871  return v_GetTimeStep(inarray);
872  }
873 
874  /**
875  * @brief Return the timestep to be used for the next step in the
876  * time-marching loop.
877  *
878  * @see UnsteadySystem::GetTimeStep
879  */
881  const Array<OneD, const Array<OneD, NekDouble> > &inarray)
882  {
883  ASSERTL0(false, "Not defined for this class");
884  return 0.0;
885  }
886 
888  {
889  return false;
890  }
891 
893  {
894  return false;
895  }
896 
898  {
899  return false;
900  }
901  }
902 }
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:135
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.
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)
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.
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.
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
std::vector< std::pair< std::string, FilterParams > > FilterMap
Definition: SessionReader.h:66
Runge-Kutta multi-stage scheme 2nd order explicit.
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
Runge-Kutta multi-stage scheme 2nd order explicit.
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)
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:1016
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