12#include "gadgetconfig.h"
14#include <gsl/gsl_linalg.h>
16#include "../data/sph_particle_data.h"
21#ifdef IMPROVED_VELOCITY_GRADIENTS
22void sph_particle_data::set_velocity_gradients(
void)
25 if(
fabs(dpos.dx_dx) > 0)
27 dvel[0][0] = dvel[0][0] / dpos.dx_dx;
35 double det = dpos.dx_dx * dpos.dy_dy - dpos.dx_dy * dpos.dx_dy;
38 double m11_inv = dpos.dy_dy / det;
39 double m12_inv = -dpos.dx_dy / det;
40 double m21_inv = -dpos.dx_dy / det;
41 double m22_inv = dpos.dx_dx / det;
43 double y11 = dvel[0][0];
44 double y21 = dvel[0][1];
45 double y12 = dvel[1][0];
46 double y22 = dvel[1][1];
48 dvel[0][0] = m11_inv * y11 + m12_inv * y21;
49 dvel[0][1] = m21_inv * y11 + m22_inv * y21;
50 dvel[1][0] = m11_inv * y12 + m12_inv * y22;
51 dvel[1][1] = m21_inv * y12 + m22_inv * y22;
52 DivVel = dvel[0][0] + dvel[1][1];
61 gsl_matrix* distance_matrix = gsl_matrix_alloc(3, 3);
62 gsl_matrix_set(distance_matrix, 0, 0, dpos.dx_dx);
63 gsl_matrix_set(distance_matrix, 0, 1, dpos.dx_dy);
64 gsl_matrix_set(distance_matrix, 0, 2, dpos.dx_dz);
65 gsl_matrix_set(distance_matrix, 1, 0, dpos.dx_dy);
66 gsl_matrix_set(distance_matrix, 1, 1, dpos.dy_dy);
67 gsl_matrix_set(distance_matrix, 1, 2, dpos.dy_dz);
68 gsl_matrix_set(distance_matrix, 2, 0, dpos.dx_dz);
69 gsl_matrix_set(distance_matrix, 2, 1, dpos.dy_dz);
70 gsl_matrix_set(distance_matrix, 2, 2, dpos.dz_dz);
73 gsl_permutation* permut = gsl_permutation_alloc(3);
74 gsl_linalg_LU_decomp(distance_matrix, permut, &sign);
76 double det = gsl_linalg_LU_det(distance_matrix, sign);
80 gsl_matrix* inv_distance_matrix = gsl_matrix_alloc(3, 3);
81 gsl_linalg_LU_invert(distance_matrix, permut, inv_distance_matrix);
84 for(
int i = 0; i < 3; i++)
86 for(
int j = 0; j < 3; j++)
88 m_inv[i][j] = gsl_matrix_get(inv_distance_matrix, i, j);
93 for(
int i = 0; i < 3; i++)
95 for(
int j = 0; j < 3; j++)
101 for(
int i = 0; i < 3; i++)
103 for(
int j = 0; j < 3; j++)
106 for(
int k = 0; k < 3; k++)
108 dvel[i][j] += m_inv[j][k] * y[k][i];
113 DivVel = dvel[0][0] + dvel[1][1] + dvel[2][2];
114 Rot[0] = dvel[2][1] - dvel[1][2];
115 Rot[1] = dvel[0][2] - dvel[2][0];
116 Rot[2] = dvel[1][0] - dvel[0][1];
119 gsl_permutation_free(permut);
120 gsl_matrix_free(distance_matrix);
121 gsl_matrix_free(inv_distance_matrix);
136#ifdef TIMEDEP_ART_VISC
137void sph_particle_data::set_viscosity_coefficient(
double dt)
139 double dDivVel_dt = dt > 0 ? (
DivVel - DivVelOld) / (dt) : 0;
141 double shockIndicator = -dDivVel_dt > 0 ? -dDivVel_dt : 0;
143 double hsml2_p = hsml_p * hsml_p;
145 double alpha_tar = (hsml2_p * shockIndicator) / (hsml2_p * shockIndicator + csnd_p * csnd_p) *
All.
ArtBulkViscConst;
149 double CsndOverHsml2 = (csnd_p / hsml_p) * (csnd_p / hsml_p);
150 double limiter = DivVel2 / (DivVel2 + CurlVel2 + 0.00001 * CsndOverHsml2);
151#ifdef NO_SHEAR_VISCOSITY_LIMITER
155 if(Alpha < alpha_tar)
157 Alpha = alpha_tar * limiter;
163 double DecayTime = 10. * hsml_p / devayVel_p;
164 Alpha = limiter * (alpha_tar + (Alpha - alpha_tar) *
exp(-dt / DecayTime));
165 if(Alpha <
All.AlphaMin)
166 Alpha =
All.AlphaMin;
global_data_all_processes All