$extrastylesheet
00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2014 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 #include "libmesh/libmesh_config.h" 00021 #ifdef LIBMESH_HAVE_SLEPC 00022 00023 #include "libmesh/diff_system.h" 00024 #include "libmesh/eigen_time_solver.h" 00025 #include "libmesh/eigen_solver.h" 00026 #include "libmesh/sparse_matrix.h" 00027 00028 namespace libMesh 00029 { 00030 00031 // Constructor 00032 EigenTimeSolver::EigenTimeSolver (sys_type& s) 00033 : Parent(s), 00034 eigen_solver (EigenSolver<Number>::build(s.comm())), 00035 tol(pow(TOLERANCE, 5./3.)), 00036 maxits(1000), 00037 n_eigenpairs_to_compute(5), 00038 n_basis_vectors_to_use(3*n_eigenpairs_to_compute), 00039 n_converged_eigenpairs(0), 00040 n_iterations_reqd(0) 00041 { 00042 libmesh_experimental(); 00043 eigen_solver->set_eigenproblem_type(GHEP);//or GNHEP 00044 eigen_solver->set_position_of_spectrum(LARGEST_MAGNITUDE); 00045 } 00046 00047 EigenTimeSolver::~EigenTimeSolver () 00048 { 00049 } 00050 00051 void EigenTimeSolver::reinit () 00052 { 00053 // empty... 00054 } 00055 00056 void EigenTimeSolver::init () 00057 { 00058 // Add matrix "B" to _system if not already there. 00059 // The user may have already added a matrix "B" before 00060 // calling the System initialization. This would be 00061 // necessary if e.g. the System originally started life 00062 // with a different type of TimeSolver and only later 00063 // had its TimeSolver changed to an EigenTimeSolver. 00064 if (!_system.have_matrix("B")) 00065 _system.add_matrix("B"); 00066 } 00067 00068 void EigenTimeSolver::solve () 00069 { 00070 // The standard implementation is basically to call: 00071 // _diff_solver->solve(); 00072 // which internally assembles (when necessary) matrices and vectors 00073 // and calls linear solver software while also doing Newton steps (see newton_solver.C) 00074 // 00075 // The element_residual and side_residual functions below control 00076 // what happens in the interior of the element assembly loops. 00077 // We have a system reference, so it's possible to call _system.assembly() 00078 // ourselves if we want to... 00079 // 00080 // Interestingly, for the EigenSolver we don't need residuals...just Jacobians. 00081 // The Jacobian should therefore always be requested, and always return 00082 // jacobian_computed as being true. 00083 00084 // The basic plan of attack is: 00085 // .) Construct the Jacobian using _system.assembly(true,true) as we 00086 // would for a steady system. Use a flag in this class to 00087 // control behavior in element_residual and side_residual 00088 // .) Swap _system.matrix to matrix "B" (be sure to add this extra matrix during init) 00089 // .) Call _system.assembly(true,true) again, using the flag in element_residual 00090 // and side_residual to only get the mass matrix terms. 00091 // .) Send A and B to Steffen's EigenSolver interface. 00092 00093 // Assemble the spatial part (matrix A) of the operator 00094 if (!this->quiet) 00095 libMesh::out << "Assembling matrix A." << std::endl; 00096 _system.matrix = &( _system.get_matrix ("System Matrix") ); 00097 this->now_assembling = Matrix_A; 00098 _system.assembly(true, true); 00099 //_system.matrix->print_matlab("matrix_A.m"); 00100 00101 // Point the system's matrix at B, call assembly again. 00102 if (!this->quiet) 00103 libMesh::out << "Assembling matrix B." << std::endl; 00104 _system.matrix = &( _system.get_matrix ("B") ); 00105 this->now_assembling = Matrix_B; 00106 _system.assembly(true, true); 00107 //_system.matrix->print_matlab("matrix_B.m"); 00108 00109 // Send matrices A, B to Steffen's SlepcEigenSolver interface 00110 //libmesh_here(); 00111 if (!this->quiet) 00112 libMesh::out << "Calling the EigenSolver." << std::endl; 00113 std::pair<unsigned int, unsigned int> solve_data = 00114 eigen_solver->solve_generalized (_system.get_matrix ("System Matrix"), 00115 _system.get_matrix ("B"), 00116 n_eigenpairs_to_compute, 00117 n_basis_vectors_to_use, 00118 tol, 00119 maxits); 00120 00121 this->n_converged_eigenpairs = solve_data.first; 00122 this->n_iterations_reqd = solve_data.second; 00123 } 00124 00125 00126 00127 bool EigenTimeSolver::element_residual(bool request_jacobian, 00128 DiffContext &context) 00129 { 00130 // The EigenTimeSolver always computes jacobians! 00131 libmesh_assert (request_jacobian); 00132 00133 context.elem_solution_rate_derivative = 1; 00134 context.elem_solution_derivative = 1; 00135 00136 // Assemble the operator for the spatial part. 00137 if (now_assembling == Matrix_A) 00138 { 00139 bool jacobian_computed = 00140 _system.element_time_derivative(request_jacobian, context); 00141 00142 // The user shouldn't compute a jacobian unless requested 00143 libmesh_assert(request_jacobian || !jacobian_computed); 00144 00145 bool jacobian_computed2 = 00146 _system.element_constraint(jacobian_computed, context); 00147 00148 // The user shouldn't compute a jacobian unless requested 00149 libmesh_assert (jacobian_computed || !jacobian_computed2); 00150 00151 return jacobian_computed && jacobian_computed2; 00152 00153 } 00154 00155 // Assemble the mass matrix operator 00156 else if (now_assembling == Matrix_B) 00157 { 00158 bool mass_jacobian_computed = 00159 _system.mass_residual(request_jacobian, context); 00160 00161 // Scale Jacobian by -1 to get positive matrix from new negative 00162 // mass_residual convention 00163 context.get_elem_jacobian() *= -1.0; 00164 00165 return mass_jacobian_computed; 00166 } 00167 00168 else 00169 libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling); 00170 } 00171 00172 00173 00174 bool EigenTimeSolver::side_residual(bool request_jacobian, 00175 DiffContext &context) 00176 { 00177 // The EigenTimeSolver always requests jacobians? 00178 //libmesh_assert (request_jacobian); 00179 00180 context.elem_solution_rate_derivative = 1; 00181 context.elem_solution_derivative = 1; 00182 00183 // Assemble the operator for the spatial part. 00184 if (now_assembling == Matrix_A) 00185 { 00186 bool jacobian_computed = 00187 _system.side_time_derivative(request_jacobian, context); 00188 00189 // The user shouldn't compute a jacobian unless requested 00190 libmesh_assert (request_jacobian || !jacobian_computed); 00191 00192 bool jacobian_computed2 = 00193 _system.side_constraint(jacobian_computed, context); 00194 00195 // The user shouldn't compute a jacobian unless requested 00196 libmesh_assert (jacobian_computed || !jacobian_computed2); 00197 00198 return jacobian_computed && jacobian_computed2; 00199 00200 } 00201 00202 // There is now a "side" equivalent for the mass matrix 00203 else if (now_assembling == Matrix_B) 00204 { 00205 bool mass_jacobian_computed = 00206 _system.side_mass_residual(request_jacobian, context); 00207 00208 // Scale Jacobian by -1 to get positive matrix from new negative 00209 // mass_residual convention 00210 context.get_elem_jacobian() *= -1.0; 00211 00212 return mass_jacobian_computed; 00213 } 00214 00215 else 00216 libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling); 00217 } 00218 00219 00220 00221 bool EigenTimeSolver::nonlocal_residual(bool request_jacobian, 00222 DiffContext &context) 00223 { 00224 // The EigenTimeSolver always requests jacobians? 00225 //libmesh_assert (request_jacobian); 00226 00227 // Assemble the operator for the spatial part. 00228 if (now_assembling == Matrix_A) 00229 { 00230 bool jacobian_computed = 00231 _system.nonlocal_time_derivative(request_jacobian, context); 00232 00233 // The user shouldn't compute a jacobian unless requested 00234 libmesh_assert (request_jacobian || !jacobian_computed); 00235 00236 bool jacobian_computed2 = 00237 _system.nonlocal_constraint(jacobian_computed, context); 00238 00239 // The user shouldn't compute a jacobian unless requested 00240 libmesh_assert (jacobian_computed || !jacobian_computed2); 00241 00242 return jacobian_computed && jacobian_computed2; 00243 00244 } 00245 00246 // There is now a "side" equivalent for the mass matrix 00247 else if (now_assembling == Matrix_B) 00248 { 00249 bool mass_jacobian_computed = 00250 _system.nonlocal_mass_residual(request_jacobian, context); 00251 00252 // Scale Jacobian by -1 to get positive matrix from new negative 00253 // mass_residual convention 00254 context.get_elem_jacobian() *= -1.0; 00255 00256 return mass_jacobian_computed; 00257 } 00258 00259 else 00260 libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling); 00261 } 00262 00263 } // namespace libMesh 00264 00265 #endif // LIBMESH_HAVE_SLEPC