param_elliptic_2d.C

00001 /*
00002  *   Copyright (c) 2004 Philippe Grandclement
00003  *
00004  *   This file is part of LORENE.
00005  *
00006  *   LORENE is free software; you can redistribute it and/or modify
00007  *   it under the terms of the GNU General Public License version 2
00008  *   as published by the Free Software Foundation.
00009  *
00010  *   LORENE is distributed in the hope that it will be useful,
00011  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00012  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00013  *   GNU General Public License for more details.
00014  *
00015  *   You should have received a copy of the GNU General Public License
00016  *   along with LORENE; if not, write to the Free Software
00017  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00018  *
00019  */
00020 
00021 char param_elliptic_2d_C[] = "$Header: /cvsroot/Lorene/C++/Source/Param_elliptic/param_elliptic_2d.C,v 1.1 2004/08/24 09:14:49 p_grandclement Exp $" ;
00022 
00023 /*
00024  * $Id: param_elliptic_2d.C,v 1.1 2004/08/24 09:14:49 p_grandclement Exp $
00025  * $Log: param_elliptic_2d.C,v $
00026  * Revision 1.1  2004/08/24 09:14:49  p_grandclement
00027  * Addition of some new operators, like Poisson in 2d... It now requieres the
00028  * GSL library to work.
00029  *
00030  * Also, the way a variable change is stored by a Param_elliptic is changed and
00031  * no longer uses Change_var but rather 2 Scalars. The codes using that feature
00032  * will requiere some modification. (It should concern only the ones about monopoles)
00033  *
00034  * 
00035  * $Header: /cvsroot/Lorene/C++/Source/Param_elliptic/param_elliptic_2d.C,v 1.1 2004/08/24 09:14:49 p_grandclement Exp $
00036  *
00037  */
00038 
00039 #include "headcpp.h"
00040 
00041 #include <math.h>
00042 #include <stdlib.h>
00043 
00044 #include "base_val.h" 
00045 #include "map.h"
00046 #include "ope_elementary.h"
00047 #include "param_elliptic.h"
00048 #include "change_var.h"
00049 #include "scalar.h"
00050 
00051 
00052 void Param_elliptic::set_poisson_2d(const Scalar& source, bool indic) {
00053 
00054   int dzpuis = source.get_dzpuis() ;
00055 
00056   if (type_map != MAP_AFF) {
00057     cout << "set_poisson_2d only defined for an affine mapping..." << endl ;
00058     abort() ;
00059   }
00060   else {
00061 
00062     int nz = get_mp().get_mg()->get_nzone() ;
00063     
00064     int nr ;
00065     double alpha, beta ;
00066     int m_quant, l_quant, base_r_1d ;
00067 
00068     int conte = 0 ;
00069     for (int l=0 ; l<nz ; l++) {
00070       
00071       nr = get_mp().get_mg()->get_nr(l) ;
00072       alpha = get_alpha (l) ;
00073       beta = get_beta (l) ;
00074 
00075       for (int k=0 ; k<get_mp().get_mg()->get_np(l)+1 ; k++)
00076     for (int j=0 ; j<get_mp().get_mg()->get_nt(l) ; j++) {
00077       if (operateurs[conte] != 0x0)     
00078         delete operateurs[conte] ;
00079         source.get_spectral_va().base.give_quant_numbers(l, k, j, m_quant, l_quant, base_r_1d) ;
00080         if (k!=1)
00081           if ((indic) || ((!indic) && (l_quant !=0)))
00082         operateurs[conte] = new Ope_poisson_2d (nr, base_r_1d, alpha, beta, l_quant, dzpuis) ;
00083           else
00084         operateurs[conte] = 0x0 ;
00085         else
00086           operateurs[conte] = 0x0 ;
00087         conte ++ ;
00088       }
00089     }
00090   }
00091 }
00092 
00093 void Param_elliptic::set_helmholtz_minus_2d(int zone, double masse, const Scalar& source) {
00094 
00095   int dzpuis = source.get_dzpuis() ;
00096   assert (masse > 0) ;
00097 
00098   if (type_map != MAP_AFF) {
00099     cout << "set_helmholtz_minus_2d only defined for an affine mapping..." << endl ;
00100     abort() ;
00101   }
00102   else {
00103 
00104     int nz = get_mp().get_mg()->get_nzone() ;
00105     if (zone == nz-1)  
00106       source.check_dzpuis(2) ; 
00107     int nr ;
00108     double alpha, beta ;
00109     int m_quant, l_quant, base_r_1d ;
00110 
00111     int conte = 0 ;
00112     for (int l=0 ; l<nz ; l++) {
00113       
00114       nr = get_mp().get_mg()->get_nr(l) ;
00115       alpha = get_alpha (l) ;
00116       beta = get_beta (l) ;
00117 
00118       for (int k=0 ; k<get_mp().get_mg()->get_np(l)+1 ; k++)
00119     for (int j=0 ; j<get_mp().get_mg()->get_nt(l) ; j++) {
00120       if (l==zone) {
00121         if (operateurs[conte] != 0x0) {     
00122           delete operateurs[conte] ;
00123           source.get_spectral_va().base.give_quant_numbers 
00124         (l, k, j, m_quant, l_quant, base_r_1d) ;
00125           operateurs[conte] = new Ope_helmholtz_minus_2d (nr, base_r_1d, alpha, beta, l_quant, masse, dzpuis) ;
00126         }
00127       } 
00128       conte ++ ;
00129     }
00130     }
00131   }
00132 }
00133 

Generated on Tue Feb 7 01:35:18 2012 for LORENE by  doxygen 1.4.6