#include <assemblysystem.h>
Public Member Functions | |
AssemblySystem () | |
~AssemblySystem () | |
void | buildGlobalMatrix (Matriz &GM, DoFHandler< DIM > &dof, EqVarType &eqvar, unsigned nQPointsPerDim, AssemblyOptimization opt=NO_OPTIMIZATION, const bool checkZero=false, double tolerance=0.0) |
template<class SparsePattern > | |
void | buildSparsityPattern (SparsePattern &pattern, DoFHandler< DIM > &dof, EqVarType &eqvar, double tolerance=0.0) |
void | buildGlobalRhs (Vetor &vRhs, DoFHandler< DIM > &dof, EqVarType &eqvar, unsigned nQPointsPerDim) |
void | addNeummanCondition (const DoFHandler< DIM > &dof, Vetor &vRhs, EqVarType &eqvar, unsigned nQPointsPerDim=2, unsigned boundary_value=0) |
void | setBoundaryValue (DoFHandler< DIM > &dof, Function3D &f, MapIntDouble &map_boundary) |
Protected Member Functions | |
void | distribute_local_to_global (Matriz &GM, FullMatrix< double > &lM, const std::vector< unsigned > &local_dof_indices, const bool checkZero, double tolerance=0.0) |
Private Attributes | |
FullMatrix< double > | m_Ml |
Vector< double > | m_Bl |
Definition at line 14 of file assemblysystem.h.
AssemblySystem< EqVarType, Matriz, Vetor >::AssemblySystem | ( | ) | [inline] |
Definition at line 65 of file assemblysystem.h.
AssemblySystem< EqVarType, Matriz, Vetor >::~AssemblySystem | ( | ) | [inline] |
Definition at line 66 of file assemblysystem.h.
void AssemblySystem< EqVarType, Matriz, Vetor >::addNeummanCondition | ( | const DoFHandler< DIM > & | dof, | |
Vetor & | vRhs, | |||
EqVarType & | eqvar, | |||
unsigned | nQPointsPerDim = 2 , |
|||
unsigned | boundary_value = 0 | |||
) | [inline] |
Definition at line 294 of file assemblysystem.h.
00295 { 00296 setValuesAtElements(eqvar,false); 00297 setValuesAtFaces(eqvar,true); 00298 00299 QGauss<DIM-1> face_quadrature(nQPointsPerDim); 00300 FEFaceValues<DIM> fe_face_values(dof.get_fe(),face_quadrature, 00301 eqvar.getUpdateFlags()|update_normal_vectors|update_q_points ); 00302 const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell; 00303 unsigned n_quadrature_points = face_quadrature.n_quadrature_points; 00304 00305 00306 //Obter referencia de campos de eqvar. 00307 unsigned &qPointFace = getRefQPointFace(eqvar); 00308 std::vector<unsigned int>& local_dof_indices = getRefDofIndices(eqvar); 00309 00310 00311 eqvar.setFaceValues(&fe_face_values); 00312 eqvar.setFeValue(NULL); 00313 local_dof_indices.resize(dofs_per_cell); 00314 00315 00316 DoFHandler<DIM>::active_cell_iterator cell = dof.begin_active(), 00317 endc = dof.end(); 00318 for (; cell!=endc; ++cell) 00319 { 00320 if (cell->at_boundary()) //Celula na fronteira ? 00321 { 00322 cell->get_dof_indices(local_dof_indices); 00323 setCell(eqvar,cell); 00324 for (unsigned face=0; face < GeometryInfo<DIM>::faces_per_cell; ++face) 00325 if (cell->face(face)->boundary_indicator() == boundary_value) 00326 { 00327 fe_face_values.reinit(cell,face); 00328 for (qPointFace=0;qPointFace<n_quadrature_points;qPointFace++) 00329 { 00330 for (unsigned i=0;i<dofs_per_cell;i++) 00331 { 00332 vRhs(local_dof_indices[i])+=eqvar.localNeumannBoundary(i); 00333 } 00334 } 00335 } 00336 } 00337 } 00338 }
void AssemblySystem< EqVarType, Matriz, Vetor >::buildGlobalMatrix | ( | Matriz & | GM, | |
DoFHandler< DIM > & | dof, | |||
EqVarType & | eqvar, | |||
unsigned | nQPointsPerDim, | |||
AssemblyOptimization | opt = NO_OPTIMIZATION , |
|||
const bool | checkZero = false , |
|||
double | tolerance = 0.0 | |||
) | [inline] |
Definition at line 68 of file assemblysystem.h.
00069 { 00070 00071 setValuesAtElements(eqvar,true); 00072 setValuesAtFaces(eqvar,false); 00073 00074 unsigned i,j; 00075 QGauss<DIM> quadrature_formula(nQPointsPerDim); 00076 FEValues<DIM> fe_values (dof.get_fe(), quadrature_formula, 00077 eqvar.getUpdateFlags()); 00078 const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell; 00079 const unsigned n_quadrature_points = quadrature_formula.n_quadrature_points; 00080 00081 //Inicializar matriz global e local. 00082 GM = 0; 00083 m_Ml.reinit(dof.get_fe().dofs_per_cell,dof.get_fe().dofs_per_cell); 00084 00085 00086 //Configurar fe values da equacao variacional 00087 eqvar.setFeValue(&fe_values); 00088 00089 //Obter referencia de alguns campos de EqVarT 00090 unsigned& qPoint = getRefQPoint(eqvar); 00091 std::vector<unsigned int>& local_dof_indices = getRefDofIndices(eqvar); 00092 00093 //Alocar vetor com dofs. 00094 local_dof_indices.resize(dof.get_fe().dofs_per_cell); 00095 00096 DoFHandler<DIM>::active_cell_iterator cell = dof.begin_active(), 00097 endc = dof.end(); 00098 00099 //Temos 3 opcoes de montagem das matrizes. 00100 //EQUAL_LOCAL_SYSTEMS 00101 if (opt == EQUAL_LOCAL_SYSTEMS) 00102 { 00103 fe_values.reinit(cell); 00104 cell->get_dof_indices(local_dof_indices); 00105 setCell(eqvar,cell); 00106 eqvar.changedCellEvent(); 00107 00108 //Construir Matriz Local. 00109 m_Ml=0.0; 00110 for (qPoint=0; qPoint<n_quadrature_points; qPoint++) 00111 { 00112 for (i=0; i<dofs_per_cell; ++i) 00113 { 00114 for (j=0; j<dofs_per_cell; ++j) 00115 m_Ml(i,j) += eqvar.localEqVar(i,j); 00116 } 00117 } 00118 for (; cell!=endc; ++cell) 00119 { 00120 cell->get_dof_indices(local_dof_indices); 00121 //Distribua os nos na matriz global 00122 this->distribute_local_to_global(GM,m_Ml,local_dof_indices,checkZero,tolerance); 00123 } 00124 return; 00125 } 00126 if (opt == EQUAL_TOPOLOGY) 00127 { 00128 fe_values.reinit(cell); 00129 for (; cell!=endc; ++cell) 00130 { 00131 cell->get_dof_indices(local_dof_indices); 00132 setCell(eqvar,cell); 00133 eqvar.changedCellEvent(); 00134 //Construir Matriz Local. 00135 m_Ml=0.0; 00136 for (qPoint=0; qPoint<n_quadrature_points; qPoint++) 00137 { 00138 for (i=0; i<dofs_per_cell; ++i) 00139 { 00140 for (j=0; j<dofs_per_cell; ++j) 00141 m_Ml(i,j) += eqvar.localEqVar(i,j); 00142 } 00143 } 00144 //Distribua os nos na matriz global 00145 this->distribute_local_to_global(GM,m_Ml,local_dof_indices,checkZero,tolerance); 00146 } 00147 return; 00148 } 00149 //Ok we dont have any processing savings 00150 //Do the usual assembly 00151 00152 for (; cell!=endc; ++cell) 00153 { 00154 fe_values.reinit(cell); 00155 cell->get_dof_indices(local_dof_indices); 00156 setCell(eqvar,cell); 00157 eqvar.changedCellEvent(); 00158 //Construir Matriz Local. 00159 m_Ml=0.0; 00160 for (qPoint=0; qPoint<n_quadrature_points; qPoint++) 00161 { 00162 for (i=0; i<dofs_per_cell; ++i) 00163 { 00164 for (j=0; j<dofs_per_cell; ++j) 00165 m_Ml(i,j) += eqvar.localEqVar(i,j); 00166 } 00167 } 00168 00169 //Distribua os nos na matriz global 00170 this->distribute_local_to_global(GM,m_Ml,local_dof_indices,checkZero,tolerance); 00171 } 00172 }
void AssemblySystem< EqVarType, Matriz, Vetor >::buildGlobalRhs | ( | Vetor & | vRhs, | |
DoFHandler< DIM > & | dof, | |||
EqVarType & | eqvar, | |||
unsigned | nQPointsPerDim | |||
) | [inline] |
Definition at line 236 of file assemblysystem.h.
00237 { 00238 unsigned i; 00239 00240 setValuesAtElements(eqvar,true); 00241 setValuesAtFaces(eqvar,false); 00242 00243 QGauss<DIM> quadrature_formula(nQPointsPerDim); 00244 FEValues<DIM> fe_values (dof.get_fe(), quadrature_formula, 00245 eqvar.getUpdateFlags()); 00246 const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell; 00247 unsigned n_quadrature_points = quadrature_formula.n_quadrature_points; 00248 00249 00250 assert(vRhs.size() == dof.n_dofs()); 00251 vRhs = 0; 00252 00253 //Obter referencia dos campos 00254 std::vector<unsigned int>& local_dof_indices = getRefDofIndices(eqvar); 00255 unsigned &qPoint = getRefQPoint(eqvar); 00256 00257 00258 eqvar.setFeValue(&fe_values); 00259 m_Bl.reinit(dofs_per_cell); 00260 local_dof_indices.resize(dofs_per_cell); 00261 00262 00263 00264 DoFHandler<DIM>::active_cell_iterator cell = dof.begin_active(), 00265 endc = dof.end(); 00266 00267 for (; cell!=endc; ++cell) 00268 { 00269 fe_values.reinit(cell); 00270 cell->get_dof_indices(local_dof_indices); 00271 setCell(eqvar,cell); 00272 eqvar.changedCellRHSEvent(); 00273 m_Bl = 0; 00274 for (qPoint=0; qPoint<n_quadrature_points; ++qPoint) 00275 { 00276 eqvar.changedQPointRHSEvent(); 00277 for (i=0; i<dofs_per_cell; ++i) 00278 { 00279 m_Bl(i) += eqvar.localRhs(i); 00280 } 00281 } 00282 00283 //Distribua os nos na matriz global 00284 for (unsigned int i=0; i<dofs_per_cell; ++i) 00285 { 00286 vRhs(local_dof_indices[i]) += m_Bl(i); 00287 } 00288 } 00289 }
void AssemblySystem< EqVarType, Matriz, Vetor >::buildSparsityPattern | ( | SparsePattern & | pattern, | |
DoFHandler< DIM > & | dof, | |||
EqVarType & | eqvar, | |||
double | tolerance = 0.0 | |||
) | [inline] |
Definition at line 175 of file assemblysystem.h.
00176 { 00177 setValuesAtElements(eqvar,true); 00178 setValuesAtFaces(eqvar,false); 00179 00180 unsigned i,j; 00181 00182 QGauss<DIM> quadrature_formula(1); 00183 FEValues<DIM> fe_values (dof.get_fe(), quadrature_formula, 00184 eqvar.getUpdateFlags()); 00185 const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell; 00186 const unsigned n_quadrature_points = quadrature_formula.n_quadrature_points; 00187 00188 //Initialize local matrix. 00189 m_Ml.reinit(dof.get_fe().dofs_per_cell,dof.get_fe().dofs_per_cell); 00190 00191 00192 //Set fe values of the variational equation 00193 eqvar.setFeValue(&fe_values); 00194 00195 //Get the reference of some fields of EqVarT 00196 unsigned& qPoint = getRefQPoint(eqvar); 00197 std::vector<unsigned int>& local_dof_indices = getRefDofIndices(eqvar); 00198 00199 //Allocate the vector of dofs indices. 00200 local_dof_indices.resize(dof.get_fe().dofs_per_cell); 00201 00202 DoFHandler<DIM>::active_cell_iterator cell = dof.begin_active(), 00203 endc = dof.end(); 00204 00205 for (; cell!=endc; ++cell) 00206 { 00207 fe_values.reinit(cell); 00208 cell->get_dof_indices(local_dof_indices); 00209 setCell(eqvar,cell); 00210 eqvar.changedCellEvent(); 00211 00212 //Construir Matriz Local. 00213 m_Ml=0.0; 00214 for (qPoint=0; qPoint<n_quadrature_points; qPoint++) 00215 { 00216 for (i=0; i<dofs_per_cell; ++i) 00217 { 00218 for (j=0; j<dofs_per_cell; ++j) 00219 m_Ml(i,j) += eqvar.localEqVar(i,j); 00220 } 00221 } 00222 for (unsigned int i=0; i<dofs_per_cell; ++i) 00223 for (unsigned int j=0; j<dofs_per_cell; ++j) 00224 { 00225 if (fabs(m_Ml(i,j)) < tolerance) 00226 { 00227 //printf("Not adding position (%d,%d) for value %g\n",local_dof_indices[i],local_dof_indices[j],m_Ml(i,j)); 00228 continue; 00229 } 00230 else 00231 pattern.add (local_dof_indices[i],local_dof_indices[j]); 00232 } 00233 } 00234 }
void AssemblySystem< EqVarType, Matriz, Vetor >::distribute_local_to_global | ( | Matriz & | GM, | |
FullMatrix< double > & | lM, | |||
const std::vector< unsigned > & | local_dof_indices, | |||
const bool | checkZero, | |||
double | tolerance = 0.0 | |||
) | [inline, protected] |
Definition at line 27 of file assemblysystem.h.
00028 { 00029 const unsigned dofs_per_cell = local_dof_indices.size(); 00030 if (checkZero) 00031 { 00032 //Distribute the global matrix`s nodes. 00033 for (unsigned int i=0; i<dofs_per_cell; ++i) 00034 for (unsigned int j=0; j<dofs_per_cell; ++j) 00035 { 00036 if (fabs(lM(i,j))<tolerance) 00037 continue; 00038 else 00039 { 00040 #ifdef DEBUG 00041 if (GM.get_sparsity_pattern().exists(local_dof_indices[i],local_dof_indices[j])) 00042 { 00043 GM.add (local_dof_indices[i],local_dof_indices[j],lM(i,j)); 00044 } 00045 else 00046 printf("Trying to put value %g in invalid entry (%d,%d)\n",lM(i,j),local_dof_indices[i],local_dof_indices[j]); 00047 #else 00048 GM.add (local_dof_indices[i],local_dof_indices[j],lM(i,j)); 00049 #endif 00050 } 00051 } 00052 } 00053 else 00054 { 00055 //Distribute the local matrix`s nodes. 00056 for (unsigned int i=0; i<dofs_per_cell; ++i) 00057 for (unsigned int j=0; j<dofs_per_cell; ++j) 00058 GM.add (local_dof_indices[i],local_dof_indices[j],lM(i,j)); 00059 } 00060 00061 }
void AssemblySystem< EqVarType, Matriz, Vetor >::setBoundaryValue | ( | DoFHandler< DIM > & | dof, | |
Function3D & | f, | |||
MapIntDouble & | map_boundary | |||
) | [inline] |
Definition at line 343 of file assemblysystem.h.
00344 { 00345 assert(dof.get_fe().n_components()== f.n_components()); 00346 const FiniteElement<DIM> &fe = dof.get_fe(); 00347 std::vector<unsigned int> face_dofs (fe.dofs_per_face, 00348 DoFHandler<DIM>::invalid_dof_index); 00349 std::vector<Point<DIM> > dof_locations (face_dofs.size(), Point<DIM>()); 00350 std::vector<double> dof_values_scalar (fe.dofs_per_face); 00351 std::vector<Point<DIM-1> > unit_support_points = fe.get_unit_face_support_points(); 00352 assert(unit_support_points.size() !=0); //Tem que ser os elementos que estou habituado 00353 assert(fe.is_primitive()); 00354 Quadrature<DIM-1> aux_quad (unit_support_points); 00355 FEFaceValues<DIM> fe_values (fe, aux_quad, update_q_points); 00356 00357 DoFHandler<DIM>::active_cell_iterator cell = dof.begin_active(), 00358 endc = dof.end(); 00359 00360 if (fe.n_components() == 1) 00361 { 00362 for (; cell!=endc; ++cell) 00363 for (unsigned int face_no = 0; face_no < GeometryInfo<DIM>::faces_per_cell;++face_no) 00364 { 00365 DoFHandler<DIM>::face_iterator face = cell->face(face_no); 00366 if (face->at_boundary()) 00367 { 00368 face->get_dof_indices (face_dofs); 00369 fe_values.reinit(cell, face_no); 00370 const std::vector<Point<DIM> > &dof_locations = fe_values.get_quadrature_points (); 00371 for (unsigned i=0;i<dof_locations.size();i++) 00372 { 00373 if (f.isInDomain(dof_locations[i])) 00374 { 00375 map_boundary[face_dofs[i]] = f(dof_locations[i]); 00376 } 00377 } 00378 } 00379 } 00380 } 00381 else //Mais de uma componente, temos mais de um grau de liberdade 00382 { 00383 for (; cell!=endc; ++cell) 00384 for (unsigned int face_no = 0; face_no < GeometryInfo<DIM>::faces_per_cell;++face_no) 00385 { 00386 DoFHandler<DIM>::face_iterator face = cell->face(face_no); 00387 if (face->at_boundary()) 00388 { 00389 face->get_dof_indices (face_dofs); 00390 fe_values.reinit(cell, face_no); 00391 const std::vector<Point<DIM> > &dof_locations = fe_values.get_quadrature_points (); 00392 for (unsigned i=0;i<dof_locations.size();i++) 00393 { 00394 unsigned component = fe.face_system_to_component_index(i).first; 00395 if (f.isInDomain(dof_locations[i],component)) 00396 { 00397 map_boundary[face_dofs[i]] = f(dof_locations[i],component); 00398 } 00399 } 00400 } 00401 } 00402 } 00403 }
Vector<double> AssemblySystem< EqVarType, Matriz, Vetor >::m_Bl [private] |
Definition at line 24 of file assemblysystem.h.
FullMatrix<double> AssemblySystem< EqVarType, Matriz, Vetor >::m_Ml [private] |
Definition at line 23 of file assemblysystem.h.