[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
mhd.cpp
1 #include <cmath>
2 #include <vector>
3 
4 #include "ADTypes.hpp"
5 
6 #include "physics.h"
7 #include "mhd.h"
8 
9 
10 namespace PHiLiP {
11 namespace Physics {
12 
13 template <int dim, int nstate, typename real>
14 std::array<real,nstate> MHD<dim,nstate,real>
16  const dealii::Point<dim,real> &pos,
17  const std::array<real,nstate> &conservative_soln,
18  const real current_time,
19  const dealii::types::global_dof_index /*cell_index*/) const
20 {
21  return source_term(pos,conservative_soln,current_time);
22 }
23 
24 template <int dim, int nstate, typename real>
25 std::array<real,nstate> MHD<dim,nstate,real>
27  const dealii::Point<dim,real> &/*pos*/,
28  const std::array<real,nstate> &/*conservative_soln*/,
29  const real /*current_time*/) const
30 {
31  std::array<real,nstate> source_term;
32  for (int s=0; s<nstate; s++) {
33  source_term[s] = 0;
34  }
35 
36  return source_term;
37 }
38 
39 //incomplete
40 template <int dim, int nstate, typename real>
41 inline std::array<real,nstate> MHD<dim,nstate,real>
42 ::convert_conservative_to_primitive ( const std::array<real,nstate> &conservative_soln ) const
43 {
44  std::array<real, nstate> primitive_soln;
45 
46  real density = conservative_soln[0];
47  dealii::Tensor<1,dim,real> vel = compute_velocities (conservative_soln);
48  real pressure = compute_pressure (conservative_soln);
49 
50  primitive_soln[0] = density;
51  for (int d=0; d<dim; ++d) {
52  primitive_soln[1+d] = vel[d];
53  }
54  primitive_soln[nstate-1] = pressure;
55  return primitive_soln;
56 }
57 
58 //incomplete
59 template <int dim, int nstate, typename real>
60 inline std::array<real,nstate> MHD<dim,nstate,real>
61 ::convert_primitive_to_conservative ( const std::array<real,nstate> &primitive_soln ) const
62 {
63 
64  const real density = primitive_soln[0];
65  const dealii::Tensor<1,dim,real> velocities = extract_velocities_from_primitive(primitive_soln);
66 
67  std::array<real, nstate> conservative_soln;
68  conservative_soln[0] = density;
69  for (int d=0; d<dim; ++d) {
70  conservative_soln[1+d] = density*velocities[d];
71  }
72  conservative_soln[nstate-1] = compute_total_energy(primitive_soln);
73 
74  return conservative_soln;
75 }
76 
77 template <int dim, int nstate, typename real>
78 inline dealii::Tensor<1,dim,real> MHD<dim,nstate,real>
79 ::compute_velocities ( const std::array<real,nstate> &conservative_soln ) const
80 {
81  const real density = conservative_soln[0];
82  dealii::Tensor<1,dim,real> vel;
83  for (int d=0; d<dim; ++d) { vel[d] = conservative_soln[1+d]/density; }
84  return vel;
85 }
86 
87 template <int dim, int nstate, typename real>
88 inline real MHD<dim,nstate,real>
89 ::compute_velocity_squared ( const dealii::Tensor<1,dim,real> &velocities ) const
90 {
91  real vel2 = 0.0;
92  for (int d=0; d<dim; d++) { vel2 = vel2 + velocities[d]*velocities[d]; }
93  return vel2;
94 }
95 
96 template <int dim, int nstate, typename real>
97 inline dealii::Tensor<1,dim,real> MHD<dim,nstate,real>
98 ::extract_velocities_from_primitive ( const std::array<real,nstate> &primitive_soln ) const
99 {
100  dealii::Tensor<1,dim,real> velocities;
101  for (int d=0; d<dim; d++) { velocities[d] = primitive_soln[1+d]; }
102  return velocities;
103 }
104 
105 
106 //incomplete
107 template <int dim, int nstate, typename real>
108 inline real MHD<dim,nstate,real>
109 ::compute_total_energy ( const std::array<real,nstate> &primitive_soln ) const
110 {
111  const real density = primitive_soln[0];
112  const real pressure = primitive_soln[nstate-1];
113  const dealii::Tensor<1,dim,real> velocities = extract_velocities_from_primitive(primitive_soln);
114  const real vel2 = compute_velocity_squared(velocities);
115 
116  const real tot_energy = pressure / gamm1 + 0.5*density*vel2;
117  return tot_energy;
118 }
119 
120 //incomplete
121 template <int dim, int nstate, typename real>
122 inline real MHD<dim,nstate,real>
123 ::compute_entropy_measure ( const std::array<real,nstate> &conservative_soln ) const
124 {
125  const real density = conservative_soln[0];
126  const real pressure = compute_pressure(conservative_soln);
127  const real entropy_measure = pressure*pow(density,-gam);
128  return entropy_measure;
129 }
130 
131 //incomplete
132 template <int dim, int nstate, typename real>
133 inline real MHD<dim,nstate,real>
134 ::compute_specific_enthalpy ( const std::array<real,nstate> &conservative_soln, const real pressure ) const
135 {
136  const real density = conservative_soln[0];
137  const real total_energy = conservative_soln[nstate-1];
138  const real specific_enthalpy = (total_energy+pressure)/density;
139  return specific_enthalpy;
140 }
141 
142 
143 template <int dim, int nstate, typename real>
144 inline real MHD<dim,nstate,real>
145 ::compute_dimensional_temperature ( const std::array<real,nstate> &primitive_soln ) const
146 {
147  const real density = primitive_soln[0];
148  const real pressure = primitive_soln[nstate-1];
149  const real temperature = gam*pressure/density;
150  return temperature;
151 }
152 
153 template <int dim, int nstate, typename real>
154 inline real MHD<dim,nstate,real>
155 ::compute_temperature ( const std::array<real,nstate> &primitive_soln ) const
156 {
157  const real dimensional_temperature = compute_dimensional_temperature(primitive_soln);
158  const real temperature = dimensional_temperature ;
159  return temperature;
160 }
161 
162 template <int dim, int nstate, typename real>
163 inline real MHD<dim,nstate,real>
164 ::compute_density_from_pressure_temperature ( const real pressure, const real temperature ) const
165 {
166  const real density = gam*pressure/temperature ;
167  return density;
168 }
169 template <int dim, int nstate, typename real>
170 inline real MHD<dim,nstate,real>
171 ::compute_temperature_from_density_pressure ( const real density, const real pressure ) const
172 {
173  const real temperature = gam*pressure/density /* mach_inf_sqr*/;
174  return temperature;
175 }
176 
177 
178 template <int dim, int nstate, typename real>
179 inline real MHD<dim,nstate,real>
180 ::compute_pressure ( const std::array<real,nstate> &conservative_soln ) const
181 {
182  const real density = conservative_soln[0];
183  //std::cout << "density " << density << std::endl;
184 
185  const real tot_energy = conservative_soln[nstate-1];
186  // std::cout << "tot_energy " << tot_energy << std::endl;
187 
188  const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
189  // std::cout << "vel1 " << vel[0] << std::endl
190  // << vel[1] << std::endl
191 // << vel[2] <<std::endl;
192 
193  const real vel2 = compute_velocity_squared(vel);
194  //std::cout << "vel ^2 " << vel2 <<std::endl;
195  real pressure = gamm1*(tot_energy - 0.5*density*vel2);
196  //std::cout << "calculated pressure is" << pressure << std::endl;
197  if(pressure<0.0) {
198  std::cout<<"Cannot compute pressure..."<<std::endl;
199  std::cout<<"density "<<density<<std::endl;
200  for(int d=0;d<dim;d++) std::cout<<"vel"<<d<<" "<<vel[d]<<std::endl;
201  std::cout<<"energy "<<tot_energy<<std::endl;
202  }
203  assert(pressure>0.0);
204  //if(pressure<1e-4) pressure = 0.01;
205  return pressure;
206 }
207 
208 template <int dim, int nstate, typename real>
209 inline real MHD<dim,nstate,real>
210 ::compute_sound ( const std::array<real,nstate> &conservative_soln ) const
211 {
212  real density = conservative_soln[0];
213  //if(density<1e-4) density = 0.01;
214  if(density<0.0) {
215  std::cout<<"density"<<density<<std::endl;
216  std::abort();
217  }
218  assert(density > 0);
219  const real pressure = compute_pressure(conservative_soln);
220  //std::cout << "pressure is" << pressure << std::endl;
221  const real sound = sqrt(pressure*gam/density);
222  //std::cout << "sound is " << sound << std::endl;
223  return sound;
224 }
225 
226 template <int dim, int nstate, typename real>
227 inline real MHD<dim,nstate,real>
228 ::compute_sound ( const real density, const real pressure ) const
229 {
230  assert(density > 0);
231  const real sound = sqrt(pressure*gam/density);
232  return sound;
233 }
234 
235 template <int dim, int nstate, typename real>
236 inline real MHD<dim,nstate,real>
237 ::compute_mach_number ( const std::array<real,nstate> &conservative_soln ) const
238 {
239  const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
240  const real velocity = sqrt(compute_velocity_squared(vel));
241  const real sound = compute_sound (conservative_soln);
242  const real mach_number = velocity/sound;
243  return mach_number;
244 }
245 
246 template <int dim, int nstate, typename real>
247 inline real MHD<dim,nstate,real>
248 ::compute_magnetic_energy (const std::array<real,nstate> &conservative_soln) const
249 {
250  real magnetic_energy = 0;
251  for (int i = 1; i <= 3; ++i)
252  magnetic_energy += 1./2. * (conservative_soln[nstate - i] * conservative_soln[nstate - i] );
253  return magnetic_energy;
254 }
255 
256 
257 // Split form functions:
258 
259 template <int dim, int nstate, typename real>
260 inline real MHD<dim,nstate,real>::
261 compute_mean_density(const std::array<real,nstate> &soln_const,
262  const std::array<real,nstate> &soln_loop) const
263 {
264  return (soln_const[0] + soln_loop[0])/2.;
265 }
266 
267 template <int dim, int nstate, typename real>
268 inline real MHD<dim,nstate,real>::
269 compute_mean_pressure(const std::array<real,nstate> &soln_const,
270  const std::array<real,nstate> &soln_loop) const
271 {
272  real pressure_const = compute_pressure(soln_const);
273  real pressure_loop = compute_pressure(soln_loop);
274  return (pressure_const + pressure_loop)/2.;
275 }
276 
277 template <int dim, int nstate, typename real>
278 inline dealii::Tensor<1,dim,real> MHD<dim,nstate,real>::
279 compute_mean_velocities(const std::array<real,nstate> &soln_const,
280  const std::array<real,nstate> &soln_loop) const
281 {
282  dealii::Tensor<1,dim,real> vel_const = compute_velocities(soln_const);
283  dealii::Tensor<1,dim,real> vel_loop = compute_velocities(soln_loop);
284  //return (vel_const + vel_loop)/2.;
285  dealii::Tensor<1,dim,real> mean_vel;
286  for (int d=0; d<0; ++d) {
287  mean_vel[d] = (vel_const[d] + vel_loop[d]) * 0.5;
288  }
289  return mean_vel;
290 }
291 
292 template <int dim, int nstate, typename real>
293 inline real MHD<dim,nstate,real>::
294 compute_mean_specific_energy(const std::array<real,nstate> &soln_const,
295  const std::array<real,nstate> &soln_loop) const
296 {
297  return ((soln_const[nstate-1]/soln_const[0]) + (soln_loop[nstate-1]/soln_loop[0]))/2.;
298 }
299 
300 
301 template <int dim, int nstate, typename real>
302 std::array<dealii::Tensor<1,dim,real>,nstate> MHD<dim,nstate,real>
303 ::convective_flux (const std::array<real,nstate> &conservative_soln) const
304 {
305  std::array<dealii::Tensor<1,dim,real>,nstate> conv_flux;
306  const real density = conservative_soln[0];
307  const real pressure = compute_pressure (conservative_soln);
308  const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
309  const real specific_total_energy = conservative_soln[nstate-1]/conservative_soln[0];
310  const real specific_total_enthalpy = specific_total_energy + pressure/density;
311  const real magnetic_energy = compute_magnetic_energy(conservative_soln);
312 
313  for (int flux_dim=0; flux_dim<dim; ++flux_dim) {
314  // Density equation
315  conv_flux[0][flux_dim] = conservative_soln[1+flux_dim];
316  // Momentum equation
317  for (int velocity_dim=0; velocity_dim<dim; ++velocity_dim){
318  conv_flux[1+velocity_dim][flux_dim] = density*vel[flux_dim]*vel[velocity_dim];
319  }
320  conv_flux[1+flux_dim][flux_dim] += pressure + magnetic_energy; // Add diagonal of pressure and magnetic energy
321  // Energy equation
322  conv_flux[nstate-4][flux_dim] = density*vel[flux_dim]*specific_total_enthalpy;
323  }
324  return conv_flux;
325 }
326 
327 template <int dim, int nstate, typename real>
328 std::array<real,nstate> MHD<dim, nstate, real>
330  const std::array<real,nstate> &conservative_soln) const
331 {
332  std::cout<<"Entropy variables for MHD hasn't been done yet."<<std::endl;
333  std::abort();
334  return conservative_soln;
335 }
336 
337 template <int dim, int nstate, typename real>
338 std::array<real,nstate> MHD<dim, nstate, real>
340  const std::array<real,nstate> &entropy_var) const
341 {
342  std::cout<<"Entropy variables for MHD hasn't been done yet."<<std::endl;
343  std::abort();
344  return entropy_var;
345 }
346 
347 template <int dim, int nstate, typename real>
348 std::array<real,nstate> MHD<dim,nstate,real>
349 ::convective_normal_flux (const std::array<real,nstate> &conservative_soln, const dealii::Tensor<1,dim,real> &normal) const
350 {
351  std::array<real, nstate> conv_normal_flux;
352  const real density = conservative_soln[0];
353  const real pressure = compute_pressure (conservative_soln);
354  const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
355  //const real normal_vel = vel*normal;
356  real normal_vel = 0.0;
357  for (int d=0; d<dim; ++d) {
358  normal_vel += vel[d]*normal[d];
359  }
360  const real total_energy = conservative_soln[nstate-1];
361  const real specific_total_enthalpy = (total_energy + pressure) / density;
362 
363  const real rhoV = density*normal_vel;
364  // Density equation
365  conv_normal_flux[0] = rhoV;
366  // Momentum equation
367  for (int velocity_dim=0; velocity_dim<dim; ++velocity_dim){
368  conv_normal_flux[1+velocity_dim] = rhoV*vel[velocity_dim] + normal[velocity_dim] * pressure;
369  }
370  // Energy equation
371  conv_normal_flux[nstate-1] = rhoV*specific_total_enthalpy;
372  return conv_normal_flux;
373 }
374 
375 template <int dim, int nstate, typename real>
376 dealii::Tensor<2,nstate,real> MHD<dim,nstate,real>
378  const std::array<real,nstate> &conservative_soln,
379  const dealii::Tensor<1,dim,real> &normal) const
380 {
381  // See Blazek Appendix A.9 p. 429-430
382  const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
383  real vel_normal = 0.0;
384  for (int d=0;d<dim;d++) { vel_normal += vel[d] * normal[d]; }
385 
386  const real vel2 = compute_velocity_squared(vel);
387  const real phi = 0.5*gamm1 * vel2;
388 
389  const real density = conservative_soln[0];
390  const real tot_energy = conservative_soln[nstate-1];
391  const real E = tot_energy / density;
392  const real a1 = gam*E-phi;
393  const real a2 = gam-1.0;
394  const real a3 = gam-2.0;
395 
396  dealii::Tensor<2,nstate,real> jacobian;
397  for (int d=0; d<dim; ++d) {
398  jacobian[0][1+d] = normal[d];
399  }
400  for (int row_dim=0; row_dim<dim; ++row_dim) {
401  jacobian[1+row_dim][0] = normal[row_dim]*phi - vel[row_dim] * vel_normal;
402  for (int col_dim=0; col_dim<dim; ++col_dim){
403  if (row_dim == col_dim) {
404  jacobian[1+row_dim][1+col_dim] = vel_normal - a3*normal[row_dim]*vel[row_dim];
405  } else {
406  jacobian[1+row_dim][1+col_dim] = normal[col_dim]*vel[row_dim] - a2*normal[row_dim]*vel[col_dim];
407  }
408  }
409  jacobian[1+row_dim][nstate-1] = normal[row_dim]*a2;
410  }
411  jacobian[nstate-1][0] = vel_normal*(phi-a1);
412  for (int d=0; d<dim; ++d){
413  jacobian[nstate-1][1+d] = normal[d]*a1 - a2*vel[d]*vel_normal;
414  }
415  jacobian[nstate-1][nstate-1] = gam*vel_normal;
416 
417  return jacobian;
418 }
419 
420 template <int dim, int nstate, typename real>
421 std::array<real,nstate> MHD<dim,nstate,real>
423  const std::array<real,nstate> &conservative_soln,
424  const dealii::Tensor<1,dim,real> &normal) const
425 {
426  const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
427  std::array<real,nstate> eig;
428  real vel_dot_n = 0.0;
429  for (int d=0;d<dim;++d) { vel_dot_n += vel[d]*normal[d]; };
430  for (int i=0; i<nstate; i++) {
431  eig[i] = vel_dot_n;
432  //eig[i] = advection_speed*normal;
433 
434  //eig[i] = 1.0;
435  //eig[i] = -1.0;
436  }
437  return eig;
438 }
439 template <int dim, int nstate, typename real>
441 ::max_convective_eigenvalue (const std::array<real,nstate> &conservative_soln) const
442 {
443  //std::cout << "going to calculate max eig" << std::endl;
444  const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
445  //std::cout << "velocities calculated" << std::endl;
446 
447  const real sound = compute_sound (conservative_soln);
448  //std::cout << "sound calculated" << std::endl;
449 
450  /*const*/ real vel2 = compute_velocity_squared(vel);
451  //std::cout << "vel2 calculated" << std::endl;
452 
453  if (vel2 < 0.0001)
454  vel2 = 0.0001;
455 
456  const real max_eig = sqrt(vel2) + sound;
457  //std::cout << "max eig calculated" << std::endl;
458 
459  return max_eig;
460 }
461 
462 template <int dim, int nstate, typename real>
464 ::max_viscous_eigenvalue (const std::array<real,nstate> &/*conservative_soln*/) const
465 {
466  return 0.0;
467 }
468 
469 template <int dim, int nstate, typename real>
470 std::array<dealii::Tensor<1,dim,real>,nstate> MHD<dim,nstate,real>
472  const std::array<real,nstate> &conservative_soln,
473  const std::array<dealii::Tensor<1,dim,real>,nstate> &solution_gradient,
474  const dealii::types::global_dof_index /*cell_index*/) const
475 {
476  return dissipative_flux(conservative_soln,solution_gradient);
477 }
478 
479 template <int dim, int nstate, typename real>
480 std::array<dealii::Tensor<1,dim,real>,nstate> MHD<dim,nstate,real>
482  const std::array<real,nstate> &/*conservative_soln*/,
483  const std::array<dealii::Tensor<1,dim,real>,nstate> &/*solution_gradient*/) const
484 {
485  std::array<dealii::Tensor<1,dim,real>,nstate> diss_flux;
486  // No dissipation
487  for (int i=0; i<nstate; i++) {
488  diss_flux[i] = 0;
489  }
490  return diss_flux;
491 }
492 
493 template <int dim, int nstate, typename real>
496  const int /*boundary_type*/,
497  const dealii::Point<dim, real> &pos,
498  const dealii::Tensor<1,dim,real> &normal_int,
499  const std::array<real,nstate> &soln_int,
500  const std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_int,
501  std::array<real,nstate> &soln_bc,
502  std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_bc) const
503 {
504  std::array<real,nstate> boundary_values;
505  std::array<dealii::Tensor<1,dim,real>,nstate> boundary_gradients;
506  for (int s=0; s<nstate; s++) {
507  boundary_values[s] = this->manufactured_solution_function->value (pos, s);
508  boundary_gradients[s] = this->manufactured_solution_function->gradient (pos, s);
509  }
510 
511  for (int istate=0; istate<nstate; ++istate) {
512 
513  std::array<real,nstate> characteristic_dot_n = convective_eigenvalues(boundary_values, normal_int);
514  const bool inflow = (characteristic_dot_n[istate] <= 0.);
515 
516  if (inflow) { // Dirichlet boundary condition
517  // soln_bc[istate] = boundary_values[istate];
518  // soln_grad_bc[istate] = soln_grad_int[istate];
519 
520  soln_bc[istate] = boundary_values[istate];
521  soln_grad_bc[istate] = soln_grad_int[istate];
522 
523  } else { // Neumann boundary condition
524  // //soln_bc[istate] = soln_int[istate];
525  // //soln_bc[istate] = boundary_values[istate];
526  // soln_bc[istate] = -soln_int[istate]+2*boundary_values[istate];
527 
528  soln_bc[istate] = soln_int[istate];
529 
530  // **************************************************************************************************************
531  // Note I don't know how to properly impose the soln_grad_bc to obtain an adjoint consistent scheme
532  // Currently, Neumann boundary conditions are only imposed for the linear advection
533  // Therefore, soln_grad_bc does not affect the solution
534  // **************************************************************************************************************
535  soln_grad_bc[istate] = soln_grad_int[istate];
536  //soln_grad_bc[istate] = boundary_gradients[istate];
537  //soln_grad_bc[istate] = -soln_grad_int[istate]+2*boundary_gradients[istate];
538  }
539  }
540 }
541 
542 //template <int dim, int nstate, typename real>
543 //void MHD<dim,nstate,real>
544 //::boundary_face_values (
545 // const int boundary_type,
546 // const dealii::Point<dim, real> &pos,
547 // const dealii::Tensor<1,dim,real> &normal_int,
548 // const std::array<real,nstate> &soln_int,
549 // const std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_int,
550 // std::array<real,nstate> &soln_bc,
551 // std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_bc) const
552 //{
553 // // NEED TO PROVIDE AS INPUT **************************************
554 // const real total_inlet_pressure = pressure_inf*pow(1.0+0.5*gamm1*mach_inf_sqr, gam/gamm1);
555 // const real total_inlet_temperature = temperature_inf*pow(total_inlet_pressure/pressure_inf, gamm1/gam);
556 //
557 // if (boundary_type == 1000) {
558 // // Manufactured solution
559 // std::array<real,nstate> conservative_boundary_values;
560 // std::array<dealii::Tensor<1,dim,real>,nstate> boundary_gradients;
561 // for (int s=0; s<nstate; s++) {
562 // conservative_boundary_values[s] = this->manufactured_solution_function.value (pos, s);
563 // boundary_gradients[s] = this->manufactured_solution_function.gradient (pos, s);
564 // }
565 // std::array<real,nstate> primitive_boundary_values = convert_conservative_to_primitive(conservative_boundary_values);
566 // for (int istate=0; istate<nstate; ++istate) {
567 //
568 // std::array<real,nstate> characteristic_dot_n = convective_eigenvalues(conservative_boundary_values, normal_int);
569 // const bool inflow = (characteristic_dot_n[istate] <= 0.);
570 //
571 // if (inflow) { // Dirichlet boundary condition
572 //
573 // soln_bc[istate] = conservative_boundary_values[istate];
574 // soln_grad_bc[istate] = soln_grad_int[istate];
575 //
576 // // Only set the pressure and velocity
577 // // primitive_boundary_values[0] = soln_int[0];;
578 // // for(int d=0;d<dim;d++){
579 // // primitive_boundary_values[1+d] = soln_int[1+d]/soln_int[0];;
580 // //}
581 // conservative_boundary_values = convert_primitive_to_conservative(primitive_boundary_values);
582 // //conservative_boundary_values[nstate-1] = soln_int[nstate-1];
583 // soln_bc[istate] = conservative_boundary_values[istate];
584 //
585 // } else { // Neumann boundary condition
586 // // soln_bc[istate] = -soln_int[istate]+2*conservative_boundary_values[istate];
587 // soln_bc[istate] = soln_int[istate];
588 //
589 // // **************************************************************************************************************
590 // // Note I don't know how to properly impose the soln_grad_bc to obtain an adjoint consistent scheme
591 // // Currently, Neumann boundary conditions are only imposed for the linear advection
592 // // Therefore, soln_grad_bc does not affect the solution
593 // // **************************************************************************************************************
594 // soln_grad_bc[istate] = soln_grad_int[istate];
595 // //soln_grad_bc[istate] = boundary_gradients[istate];
596 // //soln_grad_bc[istate] = -soln_grad_int[istate]+2*boundary_gradients[istate];
597 // }
598 //
599 // // HARDCODE DIRICHLET BC
600 // soln_bc[istate] = conservative_boundary_values[istate];
601 //
602 // }
603 // } else if (boundary_type == 1001) {
604 // // No penetration,
605 // // Given by Algorithm II of the following paper
606 // // Krivodonova, L., and Berger, M.,
607 // // “High-order accurate implementation of solid wall boundary conditions in curved geometries,”
608 // // Journal of Computational Physics, vol. 211, 2006, pp. 492–512.
609 // const std::array<real,nstate> primitive_interior_values = convert_conservative_to_primitive(soln_int);
610 //
611 // // Copy density and pressure
612 // std::array<real,nstate> primitive_boundary_values;
613 // primitive_boundary_values[0] = primitive_interior_values[0];
614 // primitive_boundary_values[nstate-1] = primitive_interior_values[nstate-1];
615 //
616 // const dealii::Tensor<1,dim,real> surface_normal = -normal_int;
617 // const dealii::Tensor<1,dim,real> velocities_int = extract_velocities_from_primitive(primitive_interior_values);
618 // const dealii::Tensor<1,dim,real> velocities_bc = velocities_int - 2.0*(velocities_int*surface_normal)*surface_normal;
619 // for (int d=0; d<dim; ++d) {
620 // primitive_boundary_values[1+d] = velocities_bc[d];
621 // }
622 //
623 // soln_bc = convert_primitive_to_conservative(primitive_boundary_values);
624 //
625 // } else if (boundary_type == 1002) {
626 // // Pressure Outflow Boundary Condition (back pressure)
627 // // Carlson 2011, sec. 2.4
628 //
629 // const real back_pressure = 0.99; // Make it as an input later on
630 //
631 // const real mach_int = compute_mach_number(soln_int);
632 // const std::array<real,nstate> primitive_interior_values = convert_conservative_to_primitive(soln_int);
633 // const real pressure_int = primitive_interior_values[nstate-1];
634 //
635 // const real radicant = 1.0+0.5*gamm1*mach_inf_sqr;
636 // const real pressure_inlet = total_inlet_pressure * pow(radicant, -gam/gamm1);
637 // const real pressure_bc = (mach_int >= 1) ? pressure_int : back_pressure*pressure_inlet;
638 // const real temperature_int = compute_temperature(primitive_interior_values);
639 //
640 // // Assign primitive boundary values
641 // std::array<real,nstate> primitive_boundary_values;
642 // primitive_boundary_values[0] = compute_density_from_pressure_temperature(pressure_bc, temperature_int);
643 // for (int d=0;d<dim;d++) { primitive_boundary_values[1+d] = primitive_interior_values[1+d]; }
644 // primitive_boundary_values[nstate-1] = pressure_bc;
645 //
646 // soln_bc = convert_primitive_to_conservative(primitive_boundary_values);
647 //
648 // // Supersonic, simply extrapolate
649 // if (mach_int > 1.0) {
650 // soln_bc = soln_int;
651 // }
652 //
653 // } else if (boundary_type == 1003) {
654 // // Inflow
655 // // Carlson 2011, sec. 2.2 & sec 2.9
656 //
657 // const std::array<real,nstate> primitive_interior_values = convert_conservative_to_primitive(soln_int);
658 //
659 // const dealii::Tensor<1,dim,real> normal = -normal_int;
660 //
661 // const real density_i = primitive_interior_values[0];
662 // const dealii::Tensor<1,dim,real> velocities_i = extract_velocities_from_primitive(primitive_interior_values);
663 // const real pressure_i = primitive_interior_values[nstate-1];
664 //
665 // const real normal_vel_i = velocities_i*normal;
666 // const real sound_i = compute_sound(soln_int);
667 // //const real mach_i = std::abs(normal_vel_i)/sound_i;
668 //
669 // //const dealii::Tensor<1,dim,real> velocities_o = velocities_inf;
670 // //const real normal_vel_o = velocities_o*normal;
671 // //const real sound_o = sound_inf;
672 // //const real mach_o = mach_inf;
673 //
674 // if(mach_inf < 1.0) {
675 // //std::cout << "Subsonic inflow, mach=" << mach_i << std::endl;
676 // // Subsonic inflow, sec 2.7
677 //
678 // // Want to solve for c_b (sound_bc), to then solve for U (velocity_magnitude_bc) and M_b (mach_bc)
679 // // Eq. 37
680 // const real riemann_pos = normal_vel_i + 2.0*sound_i/gamm1;
681 // // Could evaluate enthalpy from primitive like eq.36, but easier to use the following
682 // const real specific_total_energy = soln_int[nstate-1]/density_i;
683 // const real specific_total_enthalpy = specific_total_energy + pressure_i/density_i;
684 // // Eq. 43
685 // const real a = 1.0+2.0/gamm1;
686 // const real b = -2.0*riemann_pos;
687 // const real c = 0.5*gamm1 * (riemann_pos*riemann_pos - 2.0*specific_total_enthalpy);
688 // // Eq. 42
689 // const real term1 = -0.5*b/a;
690 // const real term2= 0.5*sqrt(b*b-4.0*a*c)/a;
691 // const real sound_bc1 = term1 + term2;
692 // const real sound_bc2 = term1 - term2;
693 // // Eq. 44
694 // const real sound_bc = std::max(sound_bc1, sound_bc2);
695 // // Eq. 45
696 // //const real velocity_magnitude_bc = 2.0*sound_bc/gamm1 - riemann_pos;
697 // const real velocity_magnitude_bc = riemann_pos - 2.0*sound_bc/gamm1;
698 // const real mach_bc = velocity_magnitude_bc/sound_bc;
699 // // Eq. 46
700 // const real radicant = 1.0+0.5*gamm1*mach_bc*mach_bc;
701 // const real pressure_bc = total_inlet_pressure * pow(radicant, -gam/gamm1);
702 // const real temperature_bc = total_inlet_temperature * pow(radicant, -1.0);
703 // //std::cout << " pressure_bc " << pressure_bc << "pressure_inf" << pressure_inf << std::endl;
704 // //std::cout << " temperature_bc " << temperature_bc << "temperature_inf" << temperature_inf << std::endl;
705 // //
706 //
707 // const real density_bc = compute_density_from_pressure_temperature(pressure_bc, temperature_bc);
708 // std::array<real,nstate> primitive_boundary_values;
709 // primitive_boundary_values[0] = density_bc;
710 // for (int d=0;d<dim;d++) { primitive_boundary_values[1+d] = velocity_magnitude_bc*normal[d]; }
711 // primitive_boundary_values[nstate-1] = pressure_bc;
712 // soln_bc = convert_primitive_to_conservative(primitive_boundary_values);
713 //
714 // //std::cout << " entropy_bc " << compute_entropy_measure(soln_bc) << "entropy_inf" << entropy_inf << std::endl;
715 //
716 // } else {
717 // // Supersonic inflow, sec 2.9
718 // // Specify all quantities through
719 // // total_inlet_pressure, total_inlet_temperature, mach_inf & angle_of_attack
720 // //std::cout << "Supersonic inflow, mach=" << mach_i << std::endl;
721 // const real radicant = 1.0+0.5*gamm1*mach_inf_sqr;
722 // const real static_inlet_pressure = total_inlet_pressure * pow(radicant, -gam/gamm1);
723 // const real static_inlet_temperature = total_inlet_temperature * pow(radicant, -1.0);
724 //
725 // const real pressure_bc = static_inlet_pressure;
726 // const real temperature_bc = static_inlet_temperature;
727 // const real density_bc = compute_density_from_pressure_temperature(pressure_bc, temperature_bc);
728 // const real sound_bc = sqrt(gam * pressure_bc / density_bc);
729 // const real velocity_magnitude_bc = mach_inf * sound_bc;
730 //
731 // // Assign primitive boundary values
732 // std::array<real,nstate> primitive_boundary_values;
733 // primitive_boundary_values[0] = density_bc;
734 // for (int d=0;d<dim;d++) { primitive_boundary_values[1+d] = -velocity_magnitude_bc*normal_int[d]; } // minus since it's inflow
735 // primitive_boundary_values[nstate-1] = pressure_bc;
736 // soln_bc = convert_primitive_to_conservative(primitive_boundary_values);
737 // //std::cout << "Inlet density : " << density_bc << std::endl;
738 // //std::cout << "Inlet vel_x : " << primitive_boundary_values[1] << std::endl;
739 // //std::cout << "Inlet vel_y : " << primitive_boundary_values[2] << std::endl;
740 // //std::cout << "Inlet pressure: " << pressure_bc << std::endl;
741 // }
742 //
743 // } else if (boundary_type == 1004) {
744 // // Farfield boundary condition
745 // const real density_bc = density_inf;
746 // const real pressure_bc = 1.0/(gam*mach_inf_sqr);
747 // std::array<real,nstate> primitive_boundary_values;
748 // primitive_boundary_values[0] = density_bc;
749 // for (int d=0;d<dim;d++) { primitive_boundary_values[1+d] = velocities_inf[d]; } // minus since it's inflow
750 // primitive_boundary_values[nstate-1] = pressure_bc;
751 // soln_bc = convert_primitive_to_conservative(primitive_boundary_values);
752 // //std::cout << "Density inf " << soln_bc[0] << std::endl;
753 // //std::cout << "momxinf " << soln_bc[1] << std::endl;
754 // //std::cout << "momxinf " << soln_bc[2] << std::endl;
755 // //std::cout << "energy inf " << soln_bc[3] << std::endl;
756 // } else{
757 // std::cout << "Invalid boundary_type: " << boundary_type << std::endl;
758 // std::abort();
759 // }
760 //}
761 //
762 //template <int dim, int nstate, typename real>
763 //dealii::Vector<double> MHD<dim,nstate,real>::post_compute_derived_quantities_vector (
764 // const dealii::Vector<double> &uh,
765 // const std::vector<dealii::Tensor<1,dim> > &duh,
766 // const std::vector<dealii::Tensor<2,dim> > &dduh,
767 // const dealii::Tensor<1,dim> &normals,
768 // const dealii::Point<dim> &evaluation_points) const
769 //{
770 // std::vector<std::string> names = post_get_names ();
771 // dealii::Vector<double> computed_quantities = PhysicsBase<dim,nstate,real>::post_compute_derived_quantities_vector ( uh, duh, dduh, normals, evaluation_points);
772 // unsigned int current_data_index = computed_quantities.size() - 1;
773 // computed_quantities.grow_or_shrink(names.size());
774 // if constexpr (std::is_same<real,double>::value) {
775 //
776 // std::array<double, nstate> conservative_soln;
777 // for (unsigned int s=0; s<nstate; ++s) {
778 // conservative_soln[s] = uh(s);
779 // }
780 // const std::array<double, nstate> primitive_soln = convert_conservative_to_primitive(conservative_soln);
781 //
782 // // Density
783 // computed_quantities(++current_data_index) = primitive_soln[0];
784 // // Velocities
785 // for (unsigned int d=0; d<dim; ++d) {
786 // computed_quantities(++current_data_index) = primitive_soln[1+d];
787 // }
788 // // Momentum
789 // for (unsigned int d=0; d<dim; ++d) {
790 // computed_quantities(++current_data_index) = conservative_soln[1+d];
791 // }
792 // // Energy
793 // computed_quantities(++current_data_index) = conservative_soln[nstate-1];
794 // // Pressure
795 // computed_quantities(++current_data_index) = primitive_soln[nstate-1];
796 // // Pressure
797 // computed_quantities(++current_data_index) = compute_temperature(primitive_soln);
798 // // Entropy generation
799 // computed_quantities(++current_data_index) = compute_entropy_measure(conservative_soln) - entropy_inf;
800 // // Mach Number
801 // computed_quantities(++current_data_index) = compute_mach_number(conservative_soln);
802 //
803 // }
804 // if (computed_quantities.size()-1 != current_data_index) {
805 // std::cout << " Did not assign a value to all the data. Missing " << computed_quantities.size() - current_data_index << " variables."
806 // << " If you added a new output variable, make sure the names and DataComponentInterpretation match the above. "
807 // << std::endl;
808 // }
809 //
810 // return computed_quantities;
811 //}
812 //
813 //template <int dim, int nstate, typename real>
814 //std::vector<dealii::DataComponentInterpretation::DataComponentInterpretation> MHD<dim,nstate,real>
815 //::post_get_data_component_interpretation () const
816 //{
817 // namespace DCI = dealii::DataComponentInterpretation;
818 // std::vector<DCI::DataComponentInterpretation> interpretation = PhysicsBase<dim,nstate,real>::post_get_data_component_interpretation (); // state variables
819 // interpretation.push_back (DCI::component_is_scalar); // Density
820 // for (unsigned int d=0; d<dim; ++d) {
821 // interpretation.push_back (DCI::component_is_part_of_vector); // Velocity
822 // }
823 // for (unsigned int d=0; d<dim; ++d) {
824 // interpretation.push_back (DCI::component_is_part_of_vector); // Momentum
825 // }
826 // interpretation.push_back (DCI::component_is_scalar); // Energy
827 // interpretation.push_back (DCI::component_is_scalar); // Pressure
828 // interpretation.push_back (DCI::component_is_scalar); // Temperature
829 // interpretation.push_back (DCI::component_is_scalar); // Entropy generation
830 // interpretation.push_back (DCI::component_is_scalar); // Mach number
831 //
832 // std::vector<std::string> names = post_get_names();
833 // if (names.size() != interpretation.size()) {
834 // std::cout << "Number of DataComponentInterpretation is not the same as number of names for output file" << std::endl;
835 // }
836 // return interpretation;
837 //}
838 //
839 //
840 //template <int dim, int nstate, typename real>
841 //std::vector<std::string> MHD<dim,nstate,real> ::post_get_names () const
842 //{
843 // std::vector<std::string> names = PhysicsBase<dim,nstate,real>::post_get_names ();
844 // names.push_back ("density");
845 // for (unsigned int d=0; d<dim; ++d) {
846 // names.push_back ("velocity");
847 // }
848 // for (unsigned int d=0; d<dim; ++d) {
849 // names.push_back ("momentum");
850 // }
851 // names.push_back ("energy");
852 // names.push_back ("pressure");
853 // names.push_back ("temperature");
854 //
855 // names.push_back ("entropy_generation");
856 // names.push_back ("mach_number");
857 // return names;
858 //}
859 //
860 //template <int dim, int nstate, typename real>
861 //dealii::UpdateFlags MHD<dim,nstate,real>
862 //::post_get_needed_update_flags () const
863 //{
864 // //return update_values | update_gradients;
865 // return dealii::update_values;
866 //}
867 
868 // Instantiate explicitly
869 template class MHD < PHILIP_DIM, 8, double >;
870 template class MHD < PHILIP_DIM, 8, FadType >;
871 template class MHD < PHILIP_DIM, 8, RadType >;
872 template class MHD < PHILIP_DIM, 8, FadFadType >;
873 template class MHD < PHILIP_DIM, 8, RadFadType >;
874 
875 } // Physics namespace
876 } // PHiLiP namespace
real compute_mean_density(const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &convervative_soln2) const
Mean density given two sets of conservative solutions.
Definition: mhd.cpp:261
real max_convective_eigenvalue(const std::array< real, nstate > &soln) const
Maximum convective eigenvalue.
Definition: mhd.cpp:441
real compute_density_from_pressure_temperature(const real pressure, const real temperature) const
Given pressure and temperature, returns NON-DIMENSIONALIZED density using free-stream non-dimensional...
Definition: mhd.cpp:164
std::array< real, nstate > convert_primitive_to_conservative(const std::array< real, nstate > &primitive_soln) const
Definition: mhd.cpp:61
std::array< real, nstate > compute_conservative_variables_from_entropy_variables(const std::array< real, nstate > &entropy_var) const
Computes the conservative variables from the entropy variables.
Definition: mhd.cpp:339
real compute_temperature_from_density_pressure(const real density, const real pressure) const
Given density and pressure, returns NON-DIMENSIONALIZED temperature using free-stream non-dimensional...
Definition: mhd.cpp:171
real compute_pressure(const std::array< real, nstate > &conservative_soln) const
Evaluate pressure from conservative variables.
Definition: mhd.cpp:180
Files for the baseline physics.
Definition: ADTypes.hpp:10
real compute_magnetic_energy(const std::array< real, nstate > &conservative_soln) const
Evaluate Magnetic Energy.
Definition: mhd.cpp:248
std::array< real, nstate > convert_conservative_to_primitive(const std::array< real, nstate > &conservative_soln) const
Definition: mhd.cpp:42
dealii::Tensor< 2, nstate, real > convective_flux_directional_jacobian(const std::array< real, nstate > &conservative_soln, const dealii::Tensor< 1, dim, real > &normal) const
Convective flux Jacobian: .
Definition: mhd.cpp:377
real compute_entropy_measure(const std::array< real, nstate > &conservative_soln) const
Evaluate entropy from conservative variables.
Definition: mhd.cpp:123
real compute_velocity_squared(const dealii::Tensor< 1, dim, real > &velocities) const
Given the velocity vector , returns the dot-product .
Definition: mhd.cpp:89
dealii::Tensor< 1, dim, real > compute_velocities(const std::array< real, nstate > &conservative_soln) const
Evaluate velocities from conservative variables.
Definition: mhd.cpp:79
dealii::Tensor< 1, dim, real > extract_velocities_from_primitive(const std::array< real, nstate > &primitive_soln) const
Given primitive variables, returns velocities.
Definition: mhd.cpp:98
real compute_sound(const std::array< real, nstate > &conservative_soln) const
Evaluate speed of sound from conservative variables.
Definition: mhd.cpp:210
real max_viscous_eigenvalue(const std::array< real, nstate > &soln) const
Maximum viscous eigenvalue.
Definition: mhd.cpp:464
void boundary_face_values(const int, const dealii::Point< dim, real > &, const dealii::Tensor< 1, dim, real > &, const std::array< real, nstate > &, const std::array< dealii::Tensor< 1, dim, real >, nstate > &, std::array< real, nstate > &, std::array< dealii::Tensor< 1, dim, real >, nstate > &) const
Evaluates boundary values and gradients on the other side of the face.
Definition: mhd.cpp:495
real compute_mean_specific_energy(const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &convervative_soln2) const
Mean specific energy given two sets of conservative solutions.
Definition: mhd.cpp:294
std::array< real, nstate > convective_eigenvalues(const std::array< real, nstate > &, const dealii::Tensor< 1, dim, real > &) const
Spectral radius of convective term Jacobian is &#39;c&#39;.
Definition: mhd.cpp:422
real compute_specific_enthalpy(const std::array< real, nstate > &conservative_soln, const real pressure) const
Evaluate pressure from conservative variables.
Definition: mhd.cpp:134
std::array< dealii::Tensor< 1, dim, real >, nstate > dissipative_flux(const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::types::global_dof_index cell_index) const
Dissipative flux: 0.
Definition: mhd.cpp:471
real compute_temperature(const std::array< real, nstate > &primitive_soln) const
Given primitive variables, returns NON-DIMENSIONALIZED temperature using free-stream non-dimensionali...
Definition: mhd.cpp:155
real compute_total_energy(const std::array< real, nstate > &primitive_soln) const
Given primitive variables, returns total energy.
Definition: mhd.cpp:109
std::array< real, nstate > compute_entropy_variables(const std::array< real, nstate > &conservative_soln) const
Computes the entropy variables.
Definition: mhd.cpp:329
real compute_dimensional_temperature(const std::array< real, nstate > &primitive_soln) const
Given primitive variables, returns DIMENSIONALIZED temperature using the equation of state...
Definition: mhd.cpp:145
dealii::Tensor< 1, dim, real > compute_mean_velocities(const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &convervative_soln2) const
Mean velocities given two sets of conservative solutions.
Definition: mhd.cpp:279
Magnetohydrodynamics (MHD) equations. Derived from PhysicsBase.
Definition: mhd.h:77
std::array< real, nstate > source_term(const dealii::Point< dim, real > &pos, const std::array< real, nstate > &conservative_soln, const real current_time, const dealii::types::global_dof_index cell_index) const
Source term is zero or depends on manufactured solution.
Definition: mhd.cpp:15
real compute_mean_pressure(const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &convervative_soln2) const
Mean pressure given two sets of conservative solutions.
Definition: mhd.cpp:269
real compute_mach_number(const std::array< real, nstate > &conservative_soln) const
Given conservative variables, returns Mach number.
Definition: mhd.cpp:237
std::array< real, nstate > convective_normal_flux(const std::array< real, nstate > &conservative_soln, const dealii::Tensor< 1, dim, real > &normal) const
Convective flux: .
Definition: mhd.cpp:349
std::array< dealii::Tensor< 1, dim, real >, nstate > convective_flux(const std::array< real, nstate > &conservative_soln) const
Convective flux: .
Definition: mhd.cpp:303