00001
00002
00003 #include "nbr_spx.h"
00004 #include "utilitaires.h"
00005 #include "graphique.h"
00006 #include "math.h"
00007 #include "metric.h"
00008 #include "param.h"
00009 #include "param_elliptic.h"
00010 #include "tensor.h"
00011 #include "unites.h"
00012 #include "excision_surf.h"
00013
00014
00015
00016
00017
00018
00019 void Excision_surf::set_expa_parab(double c_theta_lap, double c_theta_fin, Scalar& expa_fin){
00020
00021
00022
00023
00024
00025
00026 if (expa.get_spectral_va().get_etat() == 0)
00027 {
00028
00029 Scalar theta (lapse.get_mp()); theta = 3.; theta.std_spectral_base();
00030
00031 set_expa() = theta;}
00032
00033
00034
00035 Scalar thetaplus = expa;
00036
00037
00038 Scalar lapse2(thetaplus.get_mp());
00039 lapse2.annule_hard();
00040 lapse2.std_spectral_base();
00041
00042 int nt = (*lapse.get_mp().get_mg()).get_nt(1);
00043 int np = (*lapse.get_mp().get_mg()).get_np(1);
00044 for (int k=0; k<np; k++)
00045 for (int j=0; j<nt; j++) {
00046
00047
00048 lapse2.set_grid_point(0,k,j,0) = lapse.val_grid_point(1,k,j,0);
00049
00050 }
00051
00052 lapse2.std_spectral_base();
00053
00054
00055
00056
00057
00058
00059
00060
00061 Scalar ff = lapse2*(c_theta_lap*thetaplus.lapang() + c_theta_fin*(thetaplus - expa_fin));
00062
00063
00064
00065 Scalar k_1 =delta_t*ff;
00066
00067
00068 Scalar theta_int = thetaplus + k_1; theta_int.std_spectral_base();
00069
00070
00071 Scalar ff_int = lapse2*(c_theta_lap*theta_int.lapang() + c_theta_fin*(theta_int - expa_fin));
00072 ff_int.set_spectral_va().ylm();
00073
00074
00075 Scalar k_2 = delta_t*ff_int;
00076
00077
00078 Scalar bound_theta = thetaplus + k_2;
00079
00080
00081
00082 set_expa()=bound_theta;
00083 set_dt_expa()= ff_int;
00084
00085 return;
00086 }
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098 void Excision_surf::set_expa_hyperb(double alpha0, double beta0, double gamma0) {
00099
00100
00101
00102
00103
00104 Scalar thetaplus = expa;
00105 thetaplus.set_spectral_va().ylm();
00106 assert (dt_expa.get_spectral_va().get_etat() != 0);
00107 Scalar d_thetaplus = dt_expa;
00108 d_thetaplus.set_spectral_va().ylm();
00109
00111
00112 Scalar lapse2(thetaplus.get_mp());
00113 lapse2.annule_hard();
00114 lapse2.std_spectral_base();
00115
00116 int nt = (*lapse.get_mp().get_mg()).get_nt(1);
00117 int np = (*lapse.get_mp().get_mg()).get_np(1);
00118 for (int k=0; k<np; k++)
00119 for (int j=0; j<nt; j++) {
00120 lapse2.set_grid_point(0,k,j,0) = lapse.val_grid_point(1,k,j,0);
00121 }
00122
00123 lapse2.std_spectral_base();
00124
00126
00128
00130
00131
00132
00133
00134
00135 Scalar k_1 = d_thetaplus;
00136 Scalar K_1 = beta0*d_thetaplus + alpha0*d_thetaplus.lapang() + gamma0*thetaplus ;
00137
00138 thetaplus = expa + 0.5*delta_t*k_1;
00139 d_thetaplus = dt_expa + 0.5*delta_t*K_1;
00140
00141
00142 Scalar k_2 = d_thetaplus;
00143 Scalar K_2 = beta0*d_thetaplus + alpha0*d_thetaplus.lapang() + gamma0*thetaplus;
00144
00145 thetaplus = expa + 0.5*delta_t*k_2;
00146 d_thetaplus = dt_expa + 0.5*delta_t*K_2;
00147
00148
00149 Scalar k_3 = d_thetaplus;
00150 Scalar K_3 = beta0*d_thetaplus + alpha0*d_thetaplus.lapang() + gamma0*thetaplus;
00151
00152 thetaplus = expa + delta_t*k_3;
00153 d_thetaplus = dt_expa + delta_t*K_3;
00154
00155
00156 Scalar k_4 = d_thetaplus;
00157 Scalar K_4 = beta0*d_thetaplus + alpha0*d_thetaplus.lapang() + gamma0*thetaplus;
00158
00159
00160
00161 thetaplus = expa + (1./6.)*delta_t*(k_1 + 2.*k_2 + 2.*k_3 + k_4);
00162 d_thetaplus = dt_expa + (1./6.)*delta_t*(K_1 + 2.*K_2 + 2.*K_3 + K_4);
00163
00164
00165 set_expa() = thetaplus;
00166 set_dt_expa() = d_thetaplus;
00167
00168
00169
00170 return;
00171
00172 }
00173