frequency_system.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 00020 // Local includes 00021 #include "libmesh/libmesh_config.h" 00022 00023 /* 00024 * Require complex arithmetic 00025 */ 00026 #if defined(LIBMESH_USE_COMPLEX_NUMBERS) 00027 00028 00029 // C++ includes 00030 #include <cstdio> // for sprintf 00031 00032 // Local includes 00033 #include "libmesh/frequency_system.h" 00034 #include "libmesh/equation_systems.h" 00035 #include "libmesh/libmesh_logging.h" 00036 #include "libmesh/linear_solver.h" 00037 #include "libmesh/numeric_vector.h" 00038 00039 namespace libMesh 00040 { 00041 00042 00043 00044 // ------------------------------------------------------------ 00045 // FrequencySystem implementation 00046 FrequencySystem::FrequencySystem (EquationSystems& es, 00047 const std::string& name, 00048 const unsigned int number) : 00049 LinearImplicitSystem (es, name, number), 00050 solve_system (NULL), 00051 _finished_set_frequencies (false), 00052 _keep_solution_duplicates (true), 00053 _finished_init (false), 00054 _finished_assemble (false) 00055 { 00056 // default value for wave speed & fluid density 00057 //_equation_systems.parameters.set<Real>("wave speed") = 340.; 00058 //_equation_systems.parameters.set<Real>("rho") = 1.225; 00059 } 00060 00061 00062 00063 FrequencySystem::~FrequencySystem () 00064 { 00065 this->clear (); 00066 00067 // the additional matrices and vectors are cleared and zero'ed in System 00068 } 00069 00070 00071 00072 00073 void FrequencySystem::clear () 00074 { 00075 LinearImplicitSystem::clear(); 00076 00077 _finished_set_frequencies = false; 00078 _keep_solution_duplicates = true; 00079 _finished_init = false; 00080 _finished_assemble = false; 00081 00082 /* 00083 * We have to distinguish between the 00084 * simple straightforward "clear()" 00085 * and the clear that also touches the 00086 * EquationSystems parameters "current frequency" etc. 00087 * Namely, when reading from file (through equation_systems_io.C 00088 * methods), the param's are read in, then the systems. 00089 * Prior to reading a system, this system gets cleared... 00090 * And there, all the previously loaded frequency parameters 00091 * would get lost... 00092 */ 00093 } 00094 00095 00096 00097 void FrequencySystem::clear_all () 00098 { 00099 this->clear (); 00100 00101 EquationSystems& es = 00102 this->get_equation_systems(); 00103 00104 // clear frequencies in the parameters section of the 00105 // EquationSystems object 00106 if (es.parameters.have_parameter<unsigned int> ("n_frequencies")) 00107 { 00108 unsigned int n_frequencies = es.parameters.get<unsigned int>("n_frequencies"); 00109 for (unsigned int n=0; n < n_frequencies; n++) 00110 es.parameters.remove(this->form_freq_param_name(n)); 00111 es.parameters.remove("current frequency"); 00112 } 00113 } 00114 00115 00116 00117 00118 void FrequencySystem::init_data () 00119 { 00120 // initialize parent data and additional solution vectors 00121 LinearImplicitSystem::init_data(); 00122 00123 // Log how long initializing the system takes 00124 START_LOG("init()", "FrequencySystem"); 00125 00126 EquationSystems& es = 00127 this->get_equation_systems(); 00128 00129 // make sure we have frequencies to solve for 00130 if (!_finished_set_frequencies) 00131 { 00132 /* 00133 * when this system was read from file, check 00134 * if this has a "n_frequencies" parameter, 00135 * and initialize us with these. 00136 */ 00137 if (es.parameters.have_parameter<unsigned int> ("n_frequencies")) 00138 { 00139 const unsigned int n_freq = 00140 es.parameters.get<unsigned int>("n_frequencies"); 00141 00142 libmesh_assert_greater (n_freq, 0); 00143 00144 _finished_set_frequencies = true; 00145 00146 this->set_current_frequency(0); 00147 } 00148 else 00149 { 00150 libMesh::err << "ERROR: Need to set frequencies before calling init(). " << std::endl; 00151 libmesh_error(); 00152 } 00153 } 00154 00155 _finished_init = true; 00156 00157 // Stop logging init() 00158 STOP_LOG("init()", "FrequencySystem"); 00159 } 00160 00161 00162 00163 void FrequencySystem::assemble () 00164 { 00165 libmesh_assert (_finished_init); 00166 00167 if (_finished_assemble) 00168 { 00169 libMesh::err << "ERROR: Matrices already assembled." << std::endl; 00170 libmesh_error (); 00171 } 00172 00173 // Log how long assemble() takes 00174 START_LOG("assemble()", "FrequencySystem"); 00175 00176 // prepare matrix with the help of the _dof_map, 00177 // fill with sparsity pattern, initialize the 00178 // additional matrices 00179 LinearImplicitSystem::assemble(); 00180 00181 //matrix.print (); 00182 //rhs.print (); 00183 00184 _finished_assemble = true; 00185 00186 // Log how long assemble() takes 00187 STOP_LOG("assemble()", "FrequencySystem"); 00188 } 00189 00190 00191 00192 void FrequencySystem::set_frequencies_by_steps (const Real base_freq, 00193 const Real freq_step, 00194 const unsigned int n_freq, 00195 const bool allocate_solution_duplicates) 00196 { 00197 this->_keep_solution_duplicates = allocate_solution_duplicates; 00198 00199 // sanity check 00200 if (_finished_set_frequencies) 00201 { 00202 libMesh::err << "ERROR: frequencies already initialized." 00203 << std::endl; 00204 libmesh_error(); 00205 } 00206 00207 EquationSystems& es = 00208 this->get_equation_systems(); 00209 00210 // store number of frequencies as parameter 00211 es.parameters.set<unsigned int>("n_frequencies") = n_freq; 00212 00213 for (unsigned int n=0; n<n_freq; n++) 00214 { 00215 // remember frequencies as parameters 00216 es.parameters.set<Real>(this->form_freq_param_name(n)) = 00217 base_freq + n * freq_step; 00218 00219 // build storage for solution vector, if wanted 00220 if (this->_keep_solution_duplicates) 00221 this->add_vector(this->form_solu_vec_name(n)); 00222 } 00223 00224 _finished_set_frequencies = true; 00225 00226 // set the current frequency 00227 this->set_current_frequency(0); 00228 } 00229 00230 00231 00232 void FrequencySystem::set_frequencies_by_range (const Real min_freq, 00233 const Real max_freq, 00234 const unsigned int n_freq, 00235 const bool allocate_solution_duplicates) 00236 { 00237 this->_keep_solution_duplicates = allocate_solution_duplicates; 00238 00239 // sanity checks 00240 libmesh_assert_greater (max_freq, min_freq); 00241 libmesh_assert_greater (n_freq, 0); 00242 00243 if (_finished_set_frequencies) 00244 { 00245 libMesh::err << "ERROR: frequencies already initialized. " << std::endl; 00246 libmesh_error(); 00247 } 00248 00249 EquationSystems& es = 00250 this->get_equation_systems(); 00251 00252 // store number of frequencies as parameter 00253 es.parameters.set<unsigned int>("n_frequencies") = n_freq; 00254 00255 // set frequencies, build solution storage 00256 for (unsigned int n=0; n<n_freq; n++) 00257 { 00258 // remember frequencies as parameters 00259 es.parameters.set<Real>(this->form_freq_param_name(n)) = 00260 min_freq + n*(max_freq-min_freq)/(n_freq-1); 00261 00262 // build storage for solution vector, if wanted 00263 if (this->_keep_solution_duplicates) 00264 System::add_vector(this->form_solu_vec_name(n)); 00265 } 00266 00267 _finished_set_frequencies = true; 00268 00269 // set the current frequency 00270 this->set_current_frequency(0); 00271 } 00272 00273 00274 00275 void FrequencySystem::set_frequencies (const std::vector<Real>& frequencies, 00276 const bool allocate_solution_duplicates) 00277 { 00278 this->_keep_solution_duplicates = allocate_solution_duplicates; 00279 00280 // sanity checks 00281 libmesh_assert(!frequencies.empty()); 00282 00283 if (_finished_set_frequencies) 00284 { 00285 libMesh::err << "ERROR: frequencies already initialized. " << std::endl; 00286 libmesh_error(); 00287 } 00288 00289 EquationSystems& es = 00290 this->get_equation_systems(); 00291 00292 // store number of frequencies as parameter 00293 es.parameters.set<unsigned int>("n_frequencies") = frequencies.size(); 00294 00295 // set frequencies, build solution storage 00296 for (unsigned int n=0; n<frequencies.size(); n++) 00297 { 00298 // remember frequencies as parameters 00299 es.parameters.set<Real>(this->form_freq_param_name(n)) = frequencies[n]; 00300 00301 // build storage for solution vector, if wanted 00302 if (this->_keep_solution_duplicates) 00303 System::add_vector(this->form_solu_vec_name(n)); 00304 } 00305 00306 _finished_set_frequencies = true; 00307 00308 // set the current frequency 00309 this->set_current_frequency(0); 00310 } 00311 00312 00313 00314 00315 unsigned int FrequencySystem::n_frequencies () const 00316 { 00317 libmesh_assert(_finished_set_frequencies); 00318 return this->get_equation_systems().parameters.get<unsigned int>("n_frequencies"); 00319 } 00320 00321 00322 00323 void FrequencySystem::solve () 00324 { 00325 libmesh_assert_greater (this->n_frequencies(), 0); 00326 00327 // Solve for all the specified frequencies 00328 this->solve (0, this->n_frequencies()-1); 00329 } 00330 00331 00332 00333 void FrequencySystem::solve (const unsigned int n_start, 00334 const unsigned int n_stop) 00335 { 00336 // Assemble the linear system, if not already done 00337 if (!_finished_assemble) 00338 this->assemble (); 00339 00340 // the user-supplied solve method _has_ to be provided by the user 00341 libmesh_assert(solve_system); 00342 00343 // existence & range checks 00344 libmesh_assert_greater (this->n_frequencies(), 0); 00345 libmesh_assert_less (n_stop, this->n_frequencies()); 00346 00347 EquationSystems& es = 00348 this->get_equation_systems(); 00349 00350 // Get the user-specified linear solver tolerance, 00351 // the user-specified maximum # of linear solver iterations, 00352 // the user-specified wave speed 00353 const Real tol = 00354 es.parameters.get<Real>("linear solver tolerance"); 00355 const unsigned int maxits = 00356 es.parameters.get<unsigned int>("linear solver maximum iterations"); 00357 00358 00359 // // return values 00360 // std::vector< std::pair<unsigned int, Real> > vec_rval; 00361 00362 // start solver loop 00363 for (unsigned int n=n_start; n<= n_stop; n++) 00364 { 00365 // set the current frequency 00366 this->set_current_frequency(n); 00367 00368 // Call the user-supplied pre-solve method 00369 START_LOG("user_pre_solve()", "FrequencySystem"); 00370 00371 this->solve_system (es, this->name()); 00372 00373 STOP_LOG("user_pre_solve()", "FrequencySystem"); 00374 00375 00376 // Solve the linear system for this specific frequency 00377 const std::pair<unsigned int, Real> rval = 00378 linear_solver->solve (*matrix, *solution, *rhs, tol, maxits); 00379 00380 _n_linear_iterations = rval.first; 00381 _final_linear_residual = rval.second; 00382 00383 vec_rval.push_back(rval); 00384 00388 if (this->_keep_solution_duplicates) 00389 this->get_vector(this->form_solu_vec_name(n)) = *solution; 00390 } 00391 00392 // sanity check 00393 //libmesh_assert_equal_to (vec_rval.size(), (n_stop-n_start+1)); 00394 00395 //return vec_rval; 00396 } 00397 00398 00399 00400 void FrequencySystem::attach_solve_function(void fptr(EquationSystems& es, 00401 const std::string& name)) 00402 { 00403 libmesh_assert(fptr); 00404 00405 solve_system = fptr; 00406 } 00407 00408 00409 00410 void FrequencySystem::set_current_frequency(unsigned int n) 00411 { 00412 libmesh_assert_less (n, n_frequencies()); 00413 00414 EquationSystems& es = 00415 this->get_equation_systems(); 00416 00417 es.parameters.set<Real>("current frequency") = 00418 es.parameters.get<Real>(this->form_freq_param_name(n)); 00419 } 00420 00421 00422 00423 std::string FrequencySystem::form_freq_param_name(const unsigned int n) const 00424 { 00425 libmesh_assert_less (n, 9999); 00426 char buf[15]; 00427 sprintf(buf, "frequency %04d", n); 00428 return (buf); 00429 } 00430 00431 00432 00433 std::string FrequencySystem::form_solu_vec_name(const unsigned int n) const 00434 { 00435 libmesh_assert_less (n, 9999); 00436 char buf[15]; 00437 sprintf(buf, "solution %04d", n); 00438 return (buf); 00439 } 00440 00441 } // namespace libMesh 00442 00443 00444 #endif // if defined(LIBMESH_USE_COMPLEX_NUMBERS)
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:47 UTC
Hosted By: