Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 
306  cpuTime = 0.0;
307  }
308 
309  // Perform any solver-specific post-integration steps
310  if (v_PostIntegrate(step))
311  {
312  break;
313  }
314 
315  // Transform data into coefficient space
316  for (i = 0; i < nvariables; ++i)
317  {
318  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
319  m_fields[m_intVariables[i]]->FwdTrans_IterPerExp(
320  fields[i],
321  m_fields[m_intVariables[i]]->UpdateCoeffs());
322  m_fields[m_intVariables[i]]->SetPhysState(false);
323  }
324 
325  // Update filters
327  for (x = m_filters.begin(); x != m_filters.end(); ++x)
328  {
329  (*x)->Update(m_fields, m_time);
330  }
331 
332  // Write out checkpoint files
333  if ((m_checksteps && step && !((step + 1) % m_checksteps)) ||
334  doCheckTime)
335  {
337  {
338  vector<bool> transformed(nfields, false);
339  for(i = 0; i < nfields; i++)
340  {
341  if (m_fields[i]->GetWaveSpace())
342  {
343  m_fields[i]->SetWaveSpace(false);
344  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
345  m_fields[i]->UpdatePhys());
346  m_fields[i]->SetPhysState(true);
347  transformed[i] = true;
348  }
349  }
350  Checkpoint_Output(nchk++);
351  for(i = 0; i < nfields; i++)
352  {
353  if (transformed[i])
354  {
355  m_fields[i]->SetWaveSpace(true);
356  m_fields[i]->HomogeneousFwdTrans(
357  m_fields[i]->GetPhys(),
358  m_fields[i]->UpdatePhys());
359  m_fields[i]->SetPhysState(false);
360  }
361  }
362  }
363  else
364  {
365  Checkpoint_Output(nchk++);
366  }
367  doCheckTime = false;
368  }
369 
370  // Step advance
371  ++step;
372  }
373 
374  // Print out summary statistics
375  if (m_session->GetComm()->GetRank() == 0)
376  {
377  if (m_cflSafetyFactor > 0.0)
378  {
379  cout << "CFL safety factor : " << m_cflSafetyFactor << endl
380  << "CFL time-step : " << m_timestep << endl;
381  }
382  cout << "Time-integration : " << intTime << "s" << endl;
383  }
384 
385  // If homogeneous, transform back into physical space if necessary.
387  {
388  for(i = 0; i < nfields; i++)
389  {
390  if (m_fields[i]->GetWaveSpace())
391  {
392  m_fields[i]->SetWaveSpace(false);
393  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
394  m_fields[i]->UpdatePhys());
395  m_fields[i]->SetPhysState(true);
396  }
397  }
398  }
399  else
400  {
401  for(i = 0; i < nvariables; ++i)
402  {
403  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
404  m_fields[m_intVariables[i]]->SetPhysState(true);
405  }
406  }
407 
408  // Finalise filters
409  for (x = m_filters.begin(); x != m_filters.end(); ++x)
410  {
411  (*x)->Finalise(m_fields, m_time);
412  }
413 
414  // Print for 1D problems
415  if(m_spacedim == 1)
416  {
417  v_AppendOutput1D(fields);
418  }
419  }
420 
421  /**
422  * @brief Sets the initial conditions.
423  */
425  {
429  }
430 
431  /**
432  * @brief Prints a summary with some information regards the
433  * time-stepping.
434  */
436  {
438  AddSummaryItem(s, "Advection",
439  m_explicitAdvection ? "explicit" : "implicit");
440  AddSummaryItem(s, "Diffusion",
441  m_explicitDiffusion ? "explicit" : "implicit");
442 
443  if (m_session->GetSolverInfo("EQTYPE")
444  == "SteadyAdvectionDiffusionReaction")
445  {
446  AddSummaryItem(s, "Reaction",
447  m_explicitReaction ? "explicit" : "implicit");
448  }
449 
450  AddSummaryItem(s, "Time Step", m_timestep);
451  AddSummaryItem(s, "No. of Steps", m_steps);
452  AddSummaryItem(s, "Checkpoints (steps)", m_checksteps);
453  AddSummaryItem(s, "Integration Type",
455  m_intScheme->GetIntegrationMethod()]);
456  }
457 
458  /**
459  * Stores the solution in a file for 1D problems only. This method has
460  * been implemented to facilitate the post-processing for 1D problems.
461  */
463  Array<OneD, Array<OneD, NekDouble> > &solution1D)
464  {
465  // Coordinates of the quadrature points in the real physical space
466  Array<OneD,NekDouble> x(GetNpoints());
467  Array<OneD,NekDouble> y(GetNpoints());
468  Array<OneD,NekDouble> z(GetNpoints());
469  m_fields[0]->GetCoords(x, y, z);
470 
471  // Print out the solution in a txt file
472  ofstream outfile;
473  outfile.open("solution1D.txt");
474  for(int i = 0; i < GetNpoints(); i++)
475  {
476  outfile << scientific << setw (17) << setprecision(16) << x[i]
477  << " " << solution1D[0][i] << endl;
478  }
479  outfile << endl << endl;
480  outfile.close();
481  }
482 
484  Array<OneD, Array<OneD, NekDouble> > &physfield,
485  Array<OneD, Array<OneD, NekDouble> > &numflux)
486  {
487  ASSERTL0(false,
488  "This function is not implemented for this equation.");
489  }
490 
492  Array<OneD, Array<OneD, NekDouble> > &physfield,
493  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
494  Array<OneD, Array<OneD, NekDouble> > &numfluxY )
495  {
496  ASSERTL0(false,
497  "This function is not implemented for this equation.");
498  }
499 
501  const Array<OneD, Array<OneD, NekDouble> > &ufield,
502  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &uflux)
503  {
504  int i, j;
505  int nTraceNumPoints = GetTraceNpoints();
506  int nvariables = m_fields.num_elements();
507  int nqvar = uflux.num_elements();
508 
509  Array<OneD, NekDouble > Fwd (nTraceNumPoints);
510  Array<OneD, NekDouble > Bwd (nTraceNumPoints);
511  Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
512  Array<OneD, NekDouble > fluxtemp(nTraceNumPoints, 0.0);
513 
514  // Get the sign of (v \cdot n), v = an arbitrary vector
515 
516  // Evaulate upwind flux:
517  // uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
518  for (j = 0; j < nqvar; ++j)
519  {
520  for (i = 0; i < nvariables ; ++i)
521  {
522  // Compute Fwd and Bwd value of ufield of i direction
523  m_fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
524 
525  // if Vn >= 0, flux = uFwd, i.e.,
526  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uFwd
527  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uFwd
528 
529  // else if Vn < 0, flux = uBwd, i.e.,
530  // edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uBwd
531  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uBwd
532 
533  m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
534  Fwd, Bwd, fluxtemp);
535 
536  // Imposing weak boundary condition with flux
537  // if Vn >= 0, uflux = uBwd at Neumann, i.e.,
538  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uBwd
539  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uBwd
540 
541  // if Vn >= 0, uflux = uFwd at Neumann, i.e.,
542  // edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uFwd
543  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uFwd
544 
545  if(m_fields[0]->GetBndCondExpansions().num_elements())
546  {
547  WeakPenaltyforScalar(i, ufield[i], fluxtemp);
548  }
549 
550  // if Vn >= 0, flux = uFwd*(tan_{\xi}^- \cdot \vec{n}),
551  // i.e,
552  // edge::eForward, uFwd \(\tan_{\xi}^Fwd \cdot \vec{n})
553  // edge::eBackward, uFwd \(\tan_{\xi}^Bwd \cdot \vec{n})
554 
555  // else if Vn < 0, flux = uBwd*(tan_{\xi}^- \cdot \vec{n}),
556  // i.e,
557  // edge::eForward, uBwd \(\tan_{\xi}^Fwd \cdot \vec{n})
558  // edge::eBackward, uBwd \(\tan_{\xi}^Bwd \cdot \vec{n})
559 
560  Vmath::Vmul(nTraceNumPoints,
561  m_traceNormals[j], 1,
562  fluxtemp, 1,
563  uflux[j][i], 1);
564  }
565  }
566  }
567 
568 
569 
571  const Array<OneD, Array<OneD, NekDouble> > &ufield,
572  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &qfield,
573  Array<OneD, Array<OneD, NekDouble> > &qflux)
574  {
575  int nTraceNumPoints = GetTraceNpoints();
576  int nvariables = m_fields.num_elements();
577  int nqvar = qfield.num_elements();
578 
579  NekDouble C11 = 1.0;
580  Array<OneD, NekDouble > Fwd(nTraceNumPoints);
581  Array<OneD, NekDouble > Bwd(nTraceNumPoints);
582  Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
583 
584  Array<OneD, NekDouble > qFwd (nTraceNumPoints);
585  Array<OneD, NekDouble > qBwd (nTraceNumPoints);
586  Array<OneD, NekDouble > qfluxtemp(nTraceNumPoints, 0.0);
587 
588  Array<OneD, NekDouble > uterm(nTraceNumPoints);
589 
590  // Evaulate upwind flux:
591  // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
592  for (int i = 0; i < nvariables; ++i)
593  {
594  qflux[i] = Array<OneD, NekDouble> (nTraceNumPoints, 0.0);
595  for (int j = 0; j < nqvar; ++j)
596  {
597  // Compute Fwd and Bwd value of ufield of jth direction
598  m_fields[i]->GetFwdBwdTracePhys(qfield[j][i],qFwd,qBwd);
599 
600  // if Vn >= 0, flux = uFwd, i.e.,
601  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick
602  // qflux = qBwd = q+
603  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick
604  // qflux = qBwd = q-
605 
606  // else if Vn < 0, flux = uBwd, i.e.,
607  // edge::eForward, if V*n<0 <=> V*n_F<0, pick
608  // qflux = qFwd = q-
609  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick
610  // qflux = qFwd = q+
611 
612  m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
613  qBwd, qFwd,
614  qfluxtemp);
615 
616  Vmath::Vmul(nTraceNumPoints,
617  m_traceNormals[j], 1,
618  qfluxtemp, 1,
619  qfluxtemp, 1);
620 
621  // Generate Stability term = - C11 ( u- - u+ )
622  m_fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
623 
624  Vmath::Vsub(nTraceNumPoints,
625  Fwd, 1, Bwd, 1,
626  uterm, 1);
627 
628  Vmath::Smul(nTraceNumPoints,
629  -1.0 * C11, uterm, 1,
630  uterm, 1);
631 
632  // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
633  Vmath::Vadd(nTraceNumPoints,
634  uterm, 1,
635  qfluxtemp, 1,
636  qfluxtemp, 1);
637 
638  // Imposing weak boundary condition with flux
639  if (m_fields[0]->GetBndCondExpansions().num_elements())
640  {
641  WeakPenaltyforVector(i, j,
642  qfield[j][i],
643  qfluxtemp,
644  C11);
645  }
646 
647  // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
648  // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
649  // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
650  Vmath::Vadd(nTraceNumPoints,
651  qfluxtemp, 1,
652  qflux[i], 1,
653  qflux[i], 1);
654  }
655  }
656  }
657 
659  {
660  if (m_session->DefinesFunction("InitialConditions"))
661  {
662  for (int i = 0; i < m_fields.num_elements(); ++i)
663  {
665 
666  vType = m_session->GetFunctionType(
667  "InitialConditions", m_session->GetVariable(i));
668 
669  if (vType == LibUtilities::eFunctionTypeFile)
670  {
671  std::string filename
672  = m_session->GetFunctionFilename(
673  "InitialConditions", m_session->GetVariable(i));
674 
675  fs::path pfilename(filename);
676 
677  // redefine path for parallel file which is in directory
678  if(fs::is_directory(pfilename))
679  {
680  fs::path metafile("Info.xml");
681  fs::path fullpath = pfilename / metafile;
682  filename = LibUtilities::PortablePath(fullpath);
683  }
684  m_fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
685 
686  // check to see if time defined
687  if (m_fieldMetaDataMap !=
689  {
691 
692  iter = m_fieldMetaDataMap.find("Time");
693  if (iter != m_fieldMetaDataMap.end())
694  {
695  time = boost::lexical_cast<NekDouble>(
696  iter->second);
697  }
698  }
699 
700  break;
701  }
702  }
703  }
704  }
705 
707  const int var,
708  const Array<OneD, const NekDouble> &physfield,
709  Array<OneD, NekDouble> &penaltyflux,
710  NekDouble time)
711  {
712  int i, e, npoints, id1, id2;
713 
714  // Number of boundary regions
715  int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
716  int Nfps, numBDEdge;
717  int nTraceNumPoints = GetTraceNpoints();
718  int cnt = 0;
719 
720  Array<OneD, NekDouble > uplus(nTraceNumPoints);
721 
722  m_fields[var]->ExtractTracePhys(physfield, uplus);
723  for (i = 0; i < nbnd; ++i)
724  {
725  // Number of boundary expansion related to that region
726  numBDEdge = m_fields[var]->
727  GetBndCondExpansions()[i]->GetExpSize();
728 
729  // Evaluate boundary values g_D or g_N from input files
731  m_session->GetFunction("InitialConditions", 0);
732 
733  npoints = m_fields[var]->
734  GetBndCondExpansions()[i]->GetNpoints();
735 
736  Array<OneD,NekDouble> BDphysics(npoints);
737  Array<OneD,NekDouble> x0(npoints,0.0);
738  Array<OneD,NekDouble> x1(npoints,0.0);
739  Array<OneD,NekDouble> x2(npoints,0.0);
740 
741  m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
742  ifunc->Evaluate(x0,x1,x2,time,BDphysics);
743 
744  // Weakly impose boundary conditions by modifying flux values
745  for (e = 0; e < numBDEdge ; ++e)
746  {
747  // Number of points on the expansion
748  Nfps = m_fields[var]->
749  GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
750 
751  id1 = m_fields[var]->
752  GetBndCondExpansions()[i]->GetPhys_Offset(e);
753 
754  id2 = m_fields[0]->GetTrace()->
755  GetPhys_Offset(m_fields[0]->GetTraceMap()->
756  GetBndCondTraceToGlobalTraceMap(cnt++));
757 
758  // For Dirichlet boundary condition: uflux = g_D
759  if (m_fields[var]->GetBndConditions()[i]->
760  GetBoundaryConditionType() == SpatialDomains::eDirichlet)
761  {
762  Vmath::Vcopy(Nfps,
763  &BDphysics[id1], 1,
764  &penaltyflux[id2], 1);
765  }
766  // For Neumann boundary condition: uflux = u+
767  else if ((m_fields[var]->GetBndConditions()[i])->
768  GetBoundaryConditionType() == SpatialDomains::eNeumann)
769  {
770  Vmath::Vcopy(Nfps,
771  &uplus[id2], 1,
772  &penaltyflux[id2], 1);
773  }
774  }
775  }
776  }
777 
778  /**
779  * Diffusion: Imposing weak boundary condition for q with flux
780  * uflux = g_D on Dirichlet boundary condition
781  * uflux = u_Fwd on Neumann boundary condition
782  */
784  const int var,
785  const int dir,
786  const Array<OneD, const NekDouble> &physfield,
787  Array<OneD, NekDouble> &penaltyflux,
788  NekDouble C11,
789  NekDouble time)
790  {
791  int i, e, npoints, id1, id2;
792  int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
793  int numBDEdge, Nfps;
794  int nTraceNumPoints = GetTraceNpoints();
795  Array<OneD, NekDouble > uterm(nTraceNumPoints);
796  Array<OneD, NekDouble > qtemp(nTraceNumPoints);
797  int cnt = 0;
798 
799  m_fields[var]->ExtractTracePhys(physfield,qtemp);
800 
801  for (i = 0; i < nbnd; ++i)
802  {
803  numBDEdge = m_fields[var]->
804  GetBndCondExpansions()[i]->GetExpSize();
805 
806  // Evaluate boundary values g_D or g_N from input files
808  m_session->GetFunction("InitialConditions", 0);
809 
810  npoints = m_fields[var]->
811  GetBndCondExpansions()[i]->GetNpoints();
812 
813  Array<OneD,NekDouble> BDphysics(npoints);
814  Array<OneD,NekDouble> x0(npoints,0.0);
815  Array<OneD,NekDouble> x1(npoints,0.0);
816  Array<OneD,NekDouble> x2(npoints,0.0);
817 
818  m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
819  ifunc->Evaluate(x0,x1,x2,time,BDphysics);
820 
821  // Weakly impose boundary conditions by modifying flux values
822  for (e = 0; e < numBDEdge ; ++e)
823  {
824  Nfps = m_fields[var]->
825  GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
826 
827  id1 = m_fields[var]->
828  GetBndCondExpansions()[i]->GetPhys_Offset(e);
829 
830  id2 = m_fields[0]->GetTrace()->
831  GetPhys_Offset(m_fields[0]->GetTraceMap()->
832  GetBndCondTraceToGlobalTraceMap(cnt++));
833 
834  // For Dirichlet boundary condition:
835  //qflux = q+ - C_11 (u+ - g_D) (nx, ny)
836  if(m_fields[var]->GetBndConditions()[i]->
837  GetBoundaryConditionType() == SpatialDomains::eDirichlet)
838  {
839  Vmath::Vmul(Nfps,
840  &m_traceNormals[dir][id2], 1,
841  &qtemp[id2], 1,
842  &penaltyflux[id2], 1);
843  }
844  // For Neumann boundary condition: qflux = g_N
845  else if((m_fields[var]->GetBndConditions()[i])->
846  GetBoundaryConditionType() == SpatialDomains::eNeumann)
847  {
848  Vmath::Vmul(Nfps,
849  &m_traceNormals[dir][id2], 1,
850  &BDphysics[id1], 1,
851  &penaltyflux[id2], 1);
852  }
853  }
854  }
855  }
856 
857  /**
858  * @brief Return the timestep to be used for the next step in the
859  * time-marching loop.
860  *
861  * This function can be overloaded to facilitate solver which utilise a
862  * CFL (or other) parameter to determine a maximum timestep under which
863  * the problem remains stable.
864  */
866  const Array<OneD, const Array<OneD, NekDouble> > &inarray)
867  {
868  return v_GetTimeStep(inarray);
869  }
870 
871  /**
872  * @brief Return the timestep to be used for the next step in the
873  * time-marching loop.
874  *
875  * @see UnsteadySystem::GetTimeStep
876  */
878  const Array<OneD, const Array<OneD, NekDouble> > &inarray)
879  {
880  ASSERTL0(false, "Not defined for this class");
881  return 0.0;
882  }
883 
885  {
886  return false;
887  }
888 
890  {
891  return false;
892  }
893  }
894 }