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