MAST
frequency_domain_linearized_conservative_fluid_elem.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 
21 // MAST includes
25 #include "fluid/flight_condition.h"
30 #include "base/nonlinear_system.h"
31 #include "mesh/fe_base.h"
32 #include "mesh/geom_elem.h"
33 #include "base/assembly_base.h"
34 
35 
36 
39  MAST::AssemblyBase& assembly,
40  const MAST::GeomElem& elem,
41  const MAST::FlightCondition& f):
42 MAST::ConservativeFluidElementBase(sys, assembly, elem, f),
43 freq(nullptr) {
44 
45 
46 }
47 
48 
49 
50 
51 
54 
55 
56 }
57 
58 
59 
60 
61 bool
63 internal_residual (bool request_jacobian,
64  ComplexVectorX& f,
65  ComplexMatrixX& jac) {
66 
67  // first get the internal residual and Jacobian from the
68  // conservative elems, and then add to it an extra dissipation
69  // to control solution discontinuities emanating from the boundary.
70 
71  //const std::vector<Real>& JxW = _fe->get_JxW();
72  //const std::vector<std::vector<Real> >& phi = _fe->get_phi();
73  const unsigned int
74  n2 = f.size();
75  //nphi = _fe->n_shape_functions();
76 
78  f_jac_x = RealMatrixX::Zero( n2, n2),
79  fm_jac_xdot = RealMatrixX::Zero( n2, n2);
80  /*mat1_n1n1 = RealMatrixX::Zero( n1, n1),
81  mat2_n1n1 = RealMatrixX::Zero( n1, n1),
82  mat3_n1n2 = RealMatrixX::Zero( n1, n2),
83  mat4_n2n2 = RealMatrixX::Zero( n2, n2),
84  AiBi_adv = RealMatrixX::Zero( n1, n2),
85  A_sens = RealMatrixX::Zero( n1, n2),
86  LS = RealMatrixX::Zero( n1, n2),
87  LS_sens = RealMatrixX::Zero( n2, n2);*/
88 
89 
91  local_jac = ComplexMatrixX::Zero( n2, n2);
92 
93 
95  //vec1_n1 = RealVectorX::Zero(n1),
96  //vec2_n1 = RealVectorX::Zero(n1),
97  local_f = RealVectorX::Zero(n2);
98  //dc = RealVectorX::Zero(dim);
99 
100  const Complex
101  iota(0., 1.);
102 
103  Real
104  omega = 0.,
105  b_V = 0.;
106 
107  (*freq)(omega);
109 
110  // df/dx. We always need the Jacobian, since it is used to calculate
111  // the residual
113  local_f,
114  f_jac_x);
115 
116  // dfm/dxdot. We always need the Jacobian, since it is used to calculate
117  // the residual
119  local_f,
120  fm_jac_xdot,
121  f_jac_x);
122 
123  // now, combine the two to return the complex Jacobian
124 
125  local_jac = (f_jac_x.cast<Complex>() * b_V + // multiply stiffness with nondim factor
126  iota * omega * fm_jac_xdot.cast<Complex>());
127 
128  if (request_jacobian)
129  jac += local_jac;
130 
131  f += local_jac * _complex_sol;
132 
133  return request_jacobian;
134 }
135 
136 
137 
138 
139 
140 bool
143  bool request_jacobian,
144  ComplexVectorX& f,
145  ComplexMatrixX& jac) {
146 
147  // first get the internal residual and Jacobian from the
148  // conservative elems, and then add to it an extra dissipation
149  // to control solution discontinuities emanating from the boundary.
150 
151  //const std::vector<Real>& JxW = _fe->get_JxW();
152  //const std::vector<std::vector<Real> >& phi = _fe->get_phi();
153  const unsigned int
154  n2 = f.size();
155  //nphi = _fe->n_shape_functions();
156 
158  f_jac_x = RealMatrixX::Zero( n2, n2),
159  fm_jac_xdot = RealMatrixX::Zero( n2, n2);
160  /*mat1_n1n1 = RealMatrixX::Zero( n1, n1),
161  mat2_n1n1 = RealMatrixX::Zero( n1, n1),
162  mat3_n1n2 = RealMatrixX::Zero( n1, n2),
163  mat4_n2n2 = RealMatrixX::Zero( n2, n2),
164  AiBi_adv = RealMatrixX::Zero( n1, n2),
165  A_sens = RealMatrixX::Zero( n1, n2),
166  LS = RealMatrixX::Zero( n1, n2),
167  LS_sens = RealMatrixX::Zero( n2, n2);*/
168 
169 
171  local_jac = ComplexMatrixX::Zero( n2, n2),
172  local_jac_sens = ComplexMatrixX::Zero( n2, n2);
173 
174 
176  //vec1_n1 = RealVectorX::Zero(n1),
177  //vec2_n1 = RealVectorX::Zero(n1),
178  local_f = RealVectorX::Zero(n2);
179  //dc = RealVectorX::Zero(dim);
180 
181  const Complex
182  iota(0., 1.);
183 
184  Real
185  omega = 0.,
186  domega = 0.,
187  b_V = 0.;
188 
189  (*freq)(omega);
190  freq->derivative(p, domega);
192 
193 
194  // df/dx. We always need the Jacobian, since it is used to calculate
195  // the residual
197  local_f,
198  f_jac_x);
199 
200  // dfm/dxdot. We always need the Jacobian, since it is used to calculate
201  // the residual
203  local_f,
204  fm_jac_xdot,
205  f_jac_x);
206 
207  // now, combine the two to return the complex Jacobian
208 
209  local_jac = (f_jac_x.cast<Complex>() * b_V + // multiply stiffness with nondim factor
210  iota * omega * fm_jac_xdot.cast<Complex>());
211 
212  local_jac_sens = (iota * domega * fm_jac_xdot.cast<Complex>());
213 
214 
215  if (request_jacobian)
216  jac += local_jac_sens;
217 
218  f += local_jac_sens * _complex_sol + local_jac * _complex_sol_sens;
219 
220  return request_jacobian;
221 }
222 
223 
224 
225 
226 bool
228 side_external_residual (bool request_jacobian,
229  ComplexVectorX& f,
230  ComplexMatrixX& jac,
231  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc) {
232 
233  // the far-field and symmetry wall, which are assumed to be stationary will
234  // only contribute to the real part of the complex Jacobian. Hence,
235  // we will use the parent class's methods for these two. The
236  // slip wall, which may be oscillating, is implemented for this element.
237 
238  const unsigned int
239  dim = _elem.dim(),
240  n1 = dim+2,
241  n2 = f.size();
242 
244  f_jac_x = RealMatrixX::Zero(n2, n2);
245 
247  local_f = RealVectorX::Zero(n2);
248 
249  Real
250  b_V = 0.;
252 
253 
254  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>> loads;
256 
257  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>>::const_iterator
258  it = loads.begin(),
259  end = loads.end();
260 
261  for ( ; it != end; it++) {
262 
263  std::vector<MAST::BoundaryConditionBase*>::const_iterator
264  bc_it = it->second.begin(),
265  bc_end = it->second.end();
266 
267  for ( ; bc_it != bc_end; bc_it++) {
268 
269  // apply all the types of loading
270  switch ((*bc_it)->type()) {
271  case MAST::SYMMETRY_WALL: {
272 
273  f_jac_x.setZero();
274  local_f.setZero();
275 
276  // We always need the Jacobian, since it is used to
277  // calculate the residual
278 
281  local_f,
282  f_jac_x,
283  it->first,
284  **bc_it);
285 
286  // multiply jacobian with the nondimensionalizing
287  // factor (V/b for flutter analysis)
288  if (request_jacobian)
289  jac += f_jac_x.cast<Complex>() * b_V;
290  f += f_jac_x.cast<Complex>() * b_V * _complex_sol;
291  }
292  break;
293 
294  case MAST::SLIP_WALL: {
295 
296  // this calculates the Jacobian and residual contribution
297  // including the nondimensionalizing factor.
298  this->slip_wall_surface_residual(request_jacobian,
299  f,
300  jac,
301  it->first,
302  **bc_it);
303  }
304  break;
305 
306  case MAST::FAR_FIELD: {
307 
308  f_jac_x.setZero();
309  local_f.setZero();
310 
311  // We always need the Jacobian, since it is used to
312  // calculate the residual
313 
316  local_f,
317  f_jac_x,
318  it->first,
319  **bc_it);
320 
321  // multiply jacobian with the nondimensionalizing
322  // factor (V/b for flutter analysis)
323  if (request_jacobian)
324  jac += f_jac_x.cast<Complex>() * b_V;
325 
326  f += f_jac_x.cast<Complex>() * b_V * _complex_sol;
327  }
328  break;
329 
330  case MAST::DIRICHLET:
331  // nothing to be done here
332  break;
333 
334  default:
335  // not implemented yet
336  libmesh_error();
337  break;
338  }
339  }
340  }
341 
342 
343  return request_jacobian;
344 }
345 
346 
347 
348 
349 
350 
351 bool
354  bool request_jacobian,
355  ComplexVectorX& f,
356  ComplexMatrixX& jac,
357  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc) {
358 
359  // the far-field and symmetry wall, which are assumed to be stationary will
360  // only contribute to the real part of the complex Jacobian. Hence,
361  // we will use the parent class's methods for these two. The
362  // slip wall, which may be oscillating, is implemented for this element.
363 
364  const unsigned int
365  dim = _elem.dim(),
366  n1 = dim+2,
367  n2 = f.size();
368 
370  f_jac_x = RealMatrixX::Zero(n2, n2);
371 
373  local_f = RealVectorX::Zero(n2);
374 
375  Real
376  b_V = 0.;
378 
379  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>> loads;
381 
382  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>>::const_iterator
383  it = loads.begin(),
384  end = loads.end();
385 
386  for ( ; it != end; it++) {
387 
388  std::vector<MAST::BoundaryConditionBase*>::const_iterator
389  bc_it = it->second.begin(),
390  bc_end = it->second.end();
391 
392  for ( ; bc_it != bc_end; bc_it++) {
393 
394  // apply all the types of loading
395  switch ((*bc_it)->type()) {
396  case MAST::SYMMETRY_WALL: {
397 
398  f_jac_x.setZero();
399  local_f.setZero();
400 
401  // We always need the Jacobian, since it is used to
402  // calculate the residual
403 
406  local_f,
407  f_jac_x,
408  it->first,
409  **bc_it);
410 
411  // multiply jacobian with the nondimensionalizing
412  // factor (V/b for flutter analysis)
413  //if (request_jacobian)
414  // jac += f_jac_x.cast<Complex>() * b_V;
415  f += f_jac_x.cast<Complex>() * b_V * _complex_sol_sens;
416  }
417  break;
418 
419  case MAST::SLIP_WALL: {
420 
421  // this calculates the Jacobian and residual contribution
422  // including the nondimensionalizing factor.
423  this->slip_wall_surface_residual_sensitivity(p, request_jacobian,
424  f,
425  jac,
426  it->first,
427  **bc_it);
428  }
429  break;
430 
431  case MAST::FAR_FIELD: {
432 
433  f_jac_x.setZero();
434  local_f.setZero();
435 
436  // We always need the Jacobian, since it is used to
437  // calculate the residual
438 
441  local_f,
442  f_jac_x,
443  it->first,
444  **bc_it);
445 
446  // multiply jacobian with the nondimensionalizing
447  // factor (V/b for flutter analysis)
448  //if (request_jacobian)
449  // jac += f_jac_x.cast<Complex>() * b_V;
450 
451  f += f_jac_x.cast<Complex>() * b_V * _complex_sol_sens;
452  }
453  break;
454 
455  case MAST::DIRICHLET:
456  // nothing to be done here
457  break;
458 
459  default:
460  // not implemented yet
461  libmesh_error();
462  break;
463  }
464  }
465  }
466 
467  return request_jacobian;
468 }
469 
470 
471 
472 
473 
474 
475 bool
477 slip_wall_surface_residual(bool request_jacobian,
478  ComplexVectorX& f,
479  ComplexMatrixX& jac,
480  const unsigned int s,
482 
483  // inviscid boundary condition without any diffusive component
484  // conditions enforced are
485  // vi ni = wi_dot (ni + dni) - ui dni (moving slip wall with deformation)
486  // tau_ij nj = 0 (because velocity gradient at wall = 0)
487  // qi ni = 0 (since heat flux occurs only on no-slip wall and far-field bc)
488 
489  // prepare the side finite element
490  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
491 
492  const std::vector<Real> &JxW = fe->get_JxW();
493  const std::vector<libMesh::Point>& normals = fe->get_normals_for_reference_coordinate();
494  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
495 
496  const unsigned int
497  dim = _elem.dim(),
498  n1 = dim+2,
499  n2 = fe->n_shape_functions()*n1;
500 
502  vec1_n1 = RealVectorX::Zero(n1),
503  uvec = RealVectorX::Zero(3),
504  dwdot_i = RealVectorX::Zero(3),
505  ni = RealVectorX::Zero(3),
506  dni = RealVectorX::Zero(3),
507  tmp = RealVectorX::Zero(6);
508 
510  Dw_i = ComplexVectorX::Zero(3),
511  Dni = ComplexVectorX::Zero(3),
512  Duvec = ComplexVectorX::Zero(3),
513  vec2_n1 = ComplexVectorX::Zero(n1),
514  vec2_n2 = ComplexVectorX::Zero(n2),
515  flux = ComplexVectorX::Zero(n1),
516  tmp_c = ComplexVectorX::Zero(6);
517 
519  mat1_n1n1 = RealMatrixX::Zero( n1, n1),
520  mat2_n1n2 = RealMatrixX::Zero( n1, n2),
521  mat3_n2n2 = RealMatrixX::Zero( n2, n2);
522 
523  libMesh::Point pt;
525 
526  // create objects to calculate the primitive solution, flux, and Jacobian
527  MAST::PrimitiveSolution primitive_sol;
529 
530  Complex
531  iota (0., 1.),
532  Dvi_ni_freq_indep = 0.,
533  Dvi_ni_freq_dep = 0.;
534 
535  Real
536  omega = 0.,
537  b_V = 0.,
538  ui_ni_steady = 0.;
539 
540 
541  // get the surface motion object from the boundary condition object
543  *vel = nullptr;
545  *n_rot = nullptr;
546 
548  *displ_perturb = nullptr;
550  *n_rot_perturb = nullptr;
551 
552  if (bc.contains("velocity"))
553  vel = &bc.get<MAST::FieldFunction<RealVectorX> >("velocity");
554 
555  if (bc.contains("normal_rotation")) {
556 
558  tmp = bc.get<MAST::FieldFunction<RealVectorX> >("normal_rotation");
559  n_rot = dynamic_cast<MAST::NormalRotationFunctionBase<RealVectorX>*>(&tmp);
560  }
561 
562  if (bc.contains("frequency_domain_displacement")) {
563 
564  displ_perturb = &bc.get<MAST::FieldFunction<ComplexVectorX> >("frequency_domain_displacement");
565  // if displ is provided then n_rot must also be provided
566  libmesh_assert(bc.contains("frequency_domain_normal_rotation"));
567 
569  tmp = bc.get<MAST::FieldFunction<ComplexVectorX> >("frequency_domain_normal_rotation");
570  n_rot_perturb = dynamic_cast<MAST::NormalRotationFunctionBase<ComplexVectorX>*>(&tmp);
571  }
572 
573 
574 
575 
576  (*freq)(omega);
578 
579 
580  for (unsigned int qp=0; qp<JxW.size(); qp++) {
581 
582  // initialize the Bmat operator for this term
583  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
584  Bmat.right_multiply(vec1_n1, _sol); // conservative sol
585  Bmat.right_multiply(vec2_n1, _complex_sol); // perturbation sol
586 
587  // initialize the primitive solution
588  primitive_sol.zero();
589  primitive_sol.init(dim,
590  vec1_n1,
593  if_viscous());
594 
595  // initialize the small-disturbance primitive sol
596  sd_primitive_sol.zero();
597  sd_primitive_sol.init(primitive_sol, vec2_n1, if_viscous());
598 
600  // Calculation of the surface velocity term. For a
601  // steady-state system,
602  // vi (ni + dni) = wdot_i (ni + dni)
603  // or vi ni = wdot_i (ni + dni) - vi dni
604  //
605  // The first order perturbation, D(.), of this is
606  // (vi + Dvi) (ni + dni + Dni) =
607  // (wdot_i + Dwdot_i) (ni + dni + Dni)
608  // or Dvi ni = Dwdot_i (ni + dni + Dni) +
609  // wdot_i (ni + dni + Dni) -
610  // vi (ni + dni + Dni) -
611  // Dvi (dni + Dni)
612  // Neglecting HOT, this becomes
613  // Dvi ni = Dwdot_i (ni + dni) +
614  // wdot_i (ni + dni + Dni) -
615  // vi (ni + dni + Dni) -
616  // Dvi dni
617  // = Dwdot_i (ni + dni) +
618  // wdot_i Dni - vi Dni - Dvi dni
619  //
621 
622  // copy the surface normal
623  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
624  ni(i_dim) = normals[qp](i_dim);
625 
626  // now check if the surface deformation is defined and
627  // needs to be applied through transpiration boundary
628  // condition
629  primitive_sol.get_uvec(uvec);
630  sd_primitive_sol.get_duvec(Duvec);
631 
632  flux.setZero();
633 
635  // contribution from the base-flow boundary condition
637  if (vel) {
638 
639  (*vel)(qpoint[qp], _time, tmp);
640  dwdot_i = tmp.topRows(3);
641  }
642 
643  if (n_rot)
644  (*n_rot)(qpoint[qp], normals[qp], _time, dni);
645 
646  ui_ni_steady = dwdot_i.dot(ni+dni) - uvec.dot(dni);
647 
648  flux += ui_ni_steady * b_V * vec2_n1; // vi_ni dcons_flux
649  flux(n1-1) += ui_ni_steady * b_V * sd_primitive_sol.dp; // vi_ni {0,0,0,0,Dp}
650 
651  if (displ_perturb) {
652 
654  // contribution from the small-disturbance boundary condition
656  (*displ_perturb)(qpoint[qp], _time, tmp_c);
657  Dw_i = tmp_c.topRows(3);
658  (*n_rot_perturb)(qpoint[qp], normals[qp], _time, Dni);
659  }
660 
661  Dvi_ni_freq_dep = Dw_i.dot(ni+dni) * iota * omega;
662  Dvi_ni_freq_indep = (dwdot_i.cast<Complex>().dot(Dni) -
663  uvec.cast<Complex>().dot(Dni) -
664  Duvec.dot(dni));
665 
666 
667  flux += (Dvi_ni_freq_indep * b_V +
668  Dvi_ni_freq_dep) * vec1_n1; // Dvi_ni cons_flux
669  flux(n1-1) += (Dvi_ni_freq_indep * b_V +
670  Dvi_ni_freq_dep) * primitive_sol.p; // Dvi_ni {0,0,0,0,p}
671 
672  // ni Dp (only for the momentun eq)
673  flux.segment(1, dim) += ((sd_primitive_sol.dp * b_V) *
674  ni.segment(0,dim).cast<Complex>());
675 
676  Bmat.vector_mult_transpose(vec2_n2, flux);
677  f += JxW[qp] * vec2_n2;
678 
679  if ( request_jacobian ) {
680 
681  // the Jacobian contribution comes only from the dp term in the
682  // residual
684  (primitive_sol,
685  ui_ni_steady,
686  normals[qp],
687  dni,
688  mat1_n1n1);
689 
690  Bmat.left_multiply(mat2_n1n2, mat1_n1n1);
691  Bmat.right_multiply_transpose(mat3_n2n2, mat2_n1n2);
692  jac += (JxW[qp] * b_V) * mat3_n2n2.cast<Complex>() ;
693  }
694  }
695 
696  return request_jacobian;
697 }
698 
699 
700 
701 
702 
703 
704 bool
707  bool request_jacobian,
708  ComplexVectorX& f,
709  ComplexMatrixX& jac,
710  const unsigned int s,
712 
713  // inviscid boundary condition without any diffusive component
714  // conditions enforced are
715  // vi ni = wi_dot (ni + dni) - ui dni (moving slip wall with deformation)
716  // tau_ij nj = 0 (because velocity gradient at wall = 0)
717  // qi ni = 0 (since heat flux occurs only on no-slip wall and far-field bc)
718 
719  // prepare the side finite element
720  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
721 
722  const std::vector<Real> &JxW = fe->get_JxW();
723  const std::vector<libMesh::Point>& normals = fe->get_normals_for_reference_coordinate();
724  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
725 
726  const unsigned int
727  dim = _elem.dim(),
728  n1 = dim+2,
729  n2 = fe->n_shape_functions()*n1;
730 
732  vec1_n1 = RealVectorX::Zero(n1),
733  uvec = RealVectorX::Zero(3),
734  ni = RealVectorX::Zero(3),
735  dni = RealVectorX::Zero(3),
736  dwdot_i = RealVectorX::Zero(3),
737  tmp = RealVectorX::Zero(6);
738 
740  vec2_n1 = ComplexVectorX::Zero(n1),
741  vec2_n2 = ComplexVectorX::Zero(n2),
742  flux = ComplexVectorX::Zero(n1),
743  tmp_c = ComplexVectorX::Zero(6),
744  Dw_i = ComplexVectorX::Zero(3),
745  Dni = ComplexVectorX::Zero(3);
746 
748  mat1_n1n1 = RealMatrixX::Zero( n1, n1),
749  mat2_n1n2 = RealMatrixX::Zero( n1, n2),
750  mat3_n2n2 = RealMatrixX::Zero( n2, n2);
751 
752  libMesh::Point pt;
754 
755  // create objects to calculate the primitive solution, flux, and Jacobian
756  MAST::PrimitiveSolution primitive_sol;
758 
759  Complex
760  iota (0., 1.),
761  Dvi_ni_freq_indep = 0.,
762  Dvi_ni_freq_dep = 0.;
763 
764  Real
765  omega = 0.,
766  domega = 0.,
767  b_V = 0.;
768 
769 
770  // get the surface motion object from the boundary condition object
772  *displ = nullptr;
774  *n_rot = nullptr;
775 
777  *displ_perturb = nullptr;
779  *n_rot_perturb = nullptr;
780 
781 
782  if (bc.contains("displacement")) {
783 
784  displ = &bc.get<MAST::FieldFunction<RealVectorX> >("displacement");
785 
786  libmesh_assert( bc.contains("normal_rotation"));
787 
789  tmp = bc.get<MAST::FieldFunction<RealVectorX> >("normal_rotation");
790  n_rot = dynamic_cast<MAST::NormalRotationFunctionBase<RealVectorX>*>(&tmp);
791  }
792 
793 
794  if (bc.contains("frequency_domain_displacement")) {
795 
796  displ_perturb = &bc.get<MAST::FieldFunction<ComplexVectorX> >("frequency_domain_displacement");
797 
798  libmesh_assert( bc.contains("frequency_domain_normal_rotation"));
799 
801  tmp = bc.get<MAST::FieldFunction<ComplexVectorX> >("frequency_domain_normal_rotation");
802  n_rot_perturb = dynamic_cast<MAST::NormalRotationFunctionBase<ComplexVectorX>*>(&tmp);
803  }
804 
805 
806  (*freq)(omega);
807  freq->derivative(p, domega);
809 
810 
811  for (unsigned int qp=0; qp<JxW.size(); qp++) {
812 
813  // initialize the Bmat operator for this term
814  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
815  Bmat.right_multiply(vec1_n1, _sol); // conservative sol
816  Bmat.right_multiply(vec2_n1, _complex_sol); // perturbation sol
817 
818  // initialize the primitive solution
819  primitive_sol.zero();
820  primitive_sol.init(dim,
821  vec1_n1,
824  if_viscous());
825 
826  // initialize the small-disturbance primitive sol
827  sd_primitive_sol.zero();
828  sd_primitive_sol.init(primitive_sol, vec2_n1, if_viscous());
829 
831  // Calculation of the surface velocity term. For a
832  // steady-state system,
833  // vi (ni + dni) = wdot_i (ni + dni)
834  // or vi ni = wdot_i (ni + dni) - vi dni
835  //
836  // The first order perturbation, D(.), of this is
837  // (vi + Dvi) (ni + dni + Dni) =
838  // (wdot_i + Dwdot_i) (ni + dni + Dni)
839  // or Dvi ni = Dwdot_i (ni + dni + Dni) +
840  // wdot_i (ni + dni + Dni) -
841  // vi (ni + dni + Dni) -
842  // Dvi (dni + Dni)
843  // Neglecting HOT, this becomes
844  // Dvi ni = Dwdot_i (ni + dni) +
845  // wdot_i (ni + dni + Dni) -
846  // vi (ni + dni + Dni) -
847  // Dvi dni
848  //
850 
851  // copy the surface normal
852  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
853  ni(i_dim) = normals[qp](i_dim);
854 
855  // now check if the surface deformation is defined and
856  // needs to be applied through transpiration boundary
857  // condition
858  primitive_sol.get_uvec(uvec);
859 
860  if (displ) {
862  // contribution from the base-flow boundary condition
864  (*displ) (qpoint[qp], _time, tmp);
865  dwdot_i = tmp.topRows(3);
866  (*n_rot)(qpoint[qp], normals[qp], _time, dni);
867  }
868 
869  if (displ_perturb) {
870 
872  // contribution from the small-disturbance boundary condition
874  (*displ_perturb)(qpoint[qp], _time, tmp_c);
875  Dw_i = tmp_c.topRows(3);
876  (*n_rot_perturb)(qpoint[qp], normals[qp], _time, Dni);
877  }
878 
879 
880  Dvi_ni_freq_dep = Dw_i.dot(ni+dni) * iota * domega;
881 
882  flux.setZero();
883  flux += Dvi_ni_freq_dep * vec1_n1; // Dvi_ni cons_flux
884  flux(n1-1) += Dvi_ni_freq_dep * primitive_sol.p; // Dvi_ni {0,0,0,0,p}
885 
886 
887  Bmat.vector_mult_transpose(vec2_n2, flux);
888  f += JxW[qp] * vec2_n2;
889 
890  if ( request_jacobian ) {
891  libmesh_assert(false); // jac sens not implemented
892 
893  }
894  }
895 
896  return request_jacobian;
897 }
898 
899 
900 
901 
void vector_mult_transpose(T &res, const T &v) const
res = v^T * [this]
This class provides the necessary functionality for spatial discretization of the conservative fluid ...
Class defines basic operations and calculation of the small disturbance primitive variables...
bool if_viscous() const
virtual bool internal_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac)
sensitivity of internal force contribution to system residual.
const MAST::GeomElem & _elem
geometric element for which the computations are performed
Definition: elem_base.h:218
void calculate_advection_flux_jacobian_for_moving_solid_wall_boundary(const MAST::PrimitiveSolution &sol, const Real ui_ni, const libMesh::Point &nvec, const RealVectorX &dnvec, RealMatrixX &mat)
bool side_external_residual(bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc)
side external force contribution to system residual
Class defines the conversion and some basic operations on primitive fluid variables used in calculati...
void init(const MAST::PrimitiveSolution &sol, const typename VectorType< ValType >::return_type &delta_sol, bool if_viscous)
virtual bool internal_residual(bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac)
internal force contribution to system residual
ComplexVectorX _complex_sol
local solution used for frequency domain analysis
Definition: elem_base.h:250
Matrix< Complex, Dynamic, 1 > ComplexVectorX
void right_multiply(T &r, const T &m) const
[R] = [this] * [M]
virtual bool velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
libMesh::Real Real
MAST::FrequencyFunction * freq
frequency function that provides the frequency for computations.
libMesh::Complex Complex
const MAST::FlightCondition * flight_condition
This defines the surface motion for use with the nonlinear fluid solver.
virtual bool symmetry_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual void derivative(const MAST::FunctionBase &f, Real &v) const
calculates the value of the function derivative and returns it in v.
void nondimensionalizing_factor(Real &v)
Matrix< Real, Dynamic, Dynamic > RealMatrixX
Matrix< Complex, Dynamic, Dynamic > ComplexMatrixX
ComplexVectorX _complex_sol_sens
local solution used for frequency domain analysis
Definition: elem_base.h:256
void get_uvec(RealVectorX &u) const
Matrix< Real, Dynamic, 1 > RealVectorX
void get_duvec(typename VectorType< ValType >::return_type &du) const
virtual std::unique_ptr< MAST::FEBase > init_side_fe(unsigned int s, bool init_grads, int extra_quadrature_order=0) const
initializes the finite element shape function and quadrature object for the side with the order of qu...
Definition: geom_elem.cpp:159
RealVectorX _sol
local solution
Definition: elem_base.h:238
virtual bool internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual
void init(const unsigned int dim, const RealVectorX &conservative_sol, const Real cp_val, const Real cv_val, bool if_viscous)
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
void _initialize_fem_interpolation_operator(const unsigned int qp, const unsigned int dim, const MAST::FEBase &fe, MAST::FEMOperatorMatrix &Bmat)
unsigned int dim() const
Definition: geom_elem.cpp:91
GasProperty gas_property
Ambient air properties.
virtual void external_side_loads_for_quadrature_elem(std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc, std::map< unsigned int, std::vector< MAST::BoundaryConditionBase * >> &loads) const
From the given list of boundary loads, this identifies the sides of the quadrature element and the lo...
Definition: geom_elem.cpp:175
void right_multiply_transpose(T &r, const T &m) const
[R] = [this]^T * [M]
virtual bool slip_wall_surface_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
sensitivity of residual of the slip wall that may be oscillating.
void left_multiply(T &r, const T &m) const
[R] = [M] * [this]
FrequencyDomainLinearizedConservativeFluidElem(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::FlightCondition &f)
virtual bool far_field_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
const ValType & get(const std::string &nm) const
returns a constant reference to the specified function
const Real & _time
time for which system is being assembled
Definition: elem_base.h:232
virtual bool side_external_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc)
sensitivity of internal force contribution to system residual.
virtual bool slip_wall_surface_residual(bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
residual of the slip wall that may be oscillating.
bool contains(const std::string &nm) const
checks if the card contains the specified property value