GCC Code Coverage Report


Directory: ./
File: DegreeOfFreedom/dof.cpp
Date: 2024-04-14 07:32:34
Exec Total Coverage
Lines: 274 333 82.3%
Branches: 211 374 56.4%

Line Branch Exec Source
1 // ______ ______ _ _ _____ ______
2 // | ____| ____| | (_)/ ____| | ____|
3 // | |__ | |__ | | _| (___ ___| |__
4 // | __| | __| | | | |\___ \ / __| __|
5 // | | | |____| |____| |____) | (__| |____
6 // |_| |______|______|_|_____/ \___|______|
7 // Finite Elements for Life Sciences and Engineering
8 //
9 // License: LGL2.1 License
10 // FELiScE default license: LICENSE in root folder
11 //
12 // Main authors: M.A. Fernandez & J.Foulon
13 //
14
15 /*!
16 \file dof.cpp
17 \authors M.A. Fernandez & J.Foulon
18 \date 07/10/2010
19 \brief Implementation of the class Dof.
20 */
21
22 // System includes
23 #include <cmath>
24 #include <unordered_set>
25
26 // External includes
27
28 // Project includes
29 #include "DegreeOfFreedom/dof.hpp"
30
31 namespace felisce
32 {
33 /*!
34 \brief Constructor of the Dof class
35 */
36 523 Dof::Dof():
37 523 m_numDof(0),
38
1/2
✓ Branch 2 taken 523 times.
✗ Branch 3 not taken.
523 m_patternIsChanged(false)
39 523 {}
40
41 /*!
42 \brief Set m_unknown, m_variable and m_pSupportDofUnknown members.
43 \param[in] unknown complete list of unkowns.
44 \param[in] variable complete list of variables.
45 \param[in] SupportDofUnknown supports dof associate to the unknowns.
46 */
47 529 void Dof::setDof(const ListUnknown& unknown, const ListVariable& variable,
48 const std::vector<SupportDofMesh*>& pSupportDofUnknown) {
49 529 m_unknown = unknown;
50 529 m_variable = variable;
51 529 m_pSupportDofUnknown = pSupportDofUnknown;
52 529 }
53
54 /*!
55 \brief Function to get global number of dof.
56 \param[in] iEle Local number of the element.
57 \param[in] iLocDof Local number of the support dof.
58 \param[in] iVar Number of the variable
59 \param[in] iComp Number of the component of the unknown
60 \param[out] numGlobalDof Global number of the support dof.
61 */
62 1051104927 void Dof::loc2glob(felInt iEle, int iLocDof, int iVar, std::size_t iComp, felInt& numGlobalDof) const {
63 1051104927 const int iUnknown = m_variable.listIdUnknownOfVariable(iVar);
64 1051104927 int iVar2 = -1;
65
66 #ifndef NDEBUG
67
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 1051104927 times.
1051104927 if ( iUnknown == -1) {
68 std::cout << "iVar: " << iVar << std::endl;
69 FEL_ERROR("This variable isn't an unknown in this problem !\n");
70 }
71
72
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 1051104927 times.
1051104927 if ((unsigned) iVar > m_variable.size() )
73 FEL_ERROR("Not Enough variable!");
74
75
1/2
✗ Branch 2 not taken.
✓ Branch 3 taken 1051104927 times.
1051104927 if ( iComp >= m_variable[iVar].numComponent()) {
76 std::ostringstream oconv;
77 oconv << "Not enough components for this variable! ";
78 oconv << "At least " << m_variable[iVar].numComponent() + 1u << " were expected but " << iComp << " were found.";
79 FEL_ERROR(oconv.str());
80 }
81 #endif
82
83 1051104927 felInt numDofExisting = 0;
84 1051104927 felInt numSupportOfDof = 0;
85 1051104927 numGlobalDof = 0;
86 // e.g. for Navier Stokes model, we have two unknowns. In m_pSupportDofUnknown, there are first the velocity dof support and then pressure one.
87 // So if iUnknown==0 (velocity), we directly compute the global dof number. Else if iUnknown==1 (pressure), we have to skip all the velocity global dof numbers (numDofExisting).
88
2/2
✓ Branch 0 taken 75487307 times.
✓ Branch 1 taken 1051104927 times.
1126592234 for ( int i = 0; i < iUnknown; i++) {
89 75487307 iVar2 = m_unknown.idVariable(i);
90 75487307 numSupportOfDof = m_pSupportDofUnknown[i]->numSupportDof();
91 75487307 numDofExisting += numSupportOfDof * m_variable[iVar2].numComponent() ;
92 }
93
94 1051104927 const auto& r_support_unknowns = m_pSupportDofUnknown[iUnknown];
95 1051104927 numGlobalDof = numDofExisting + r_support_unknowns->iSupportDof()[r_support_unknowns->iEle()[iEle] + iLocDof ] * m_variable[iVar].numComponent() + iComp;
96 1051104927 }
97
98 /*!
99 \brief Function to get global number of dof.
100 */
101 48294905 void Dof::loc2glob(felInt iEle, int numLocDof, int iVar, std::size_t numComp, std::vector<felInt>& numGlobalDof) const {
102
1/2
✓ Branch 2 taken 48294905 times.
✗ Branch 3 not taken.
48294905 std::vector<std::size_t> loc2globComp(numComp);
103
1/2
✓ Branch 2 taken 48294905 times.
✗ Branch 3 not taken.
48294905 std::vector<int> loc2globDof(numLocDof);
104
105
2/2
✓ Branch 0 taken 98723483 times.
✓ Branch 1 taken 48294905 times.
147018388 for (std::size_t iComp = 0; iComp < numComp; iComp++)
106 98723483 loc2globComp[iComp] = iComp;
107
108
2/2
✓ Branch 0 taken 185790699 times.
✓ Branch 1 taken 48294905 times.
234085604 for (int iSupport = 0; iSupport < numLocDof; iSupport++)
109 185790699 loc2globDof[iSupport] = iSupport;
110
111
1/2
✓ Branch 1 taken 48294905 times.
✗ Branch 2 not taken.
48294905 this->loc2glob(iEle, loc2globDof, iVar, loc2globComp, numGlobalDof);
112 48294905 }
113
114 /*!
115 \brief Function to get global number of dof.
116 */
117 48294905 void Dof::loc2glob(felInt iEle, std::vector<int>& iLocDof, int iVar, std::vector<std::size_t>& iComp, std::vector<felInt>& numGlobalDof) const {
118 48294905 int iUnknown = m_variable.listIdUnknownOfVariable(iVar);
119 48294905 int iVar2 = -1;
120
1/2
✓ Branch 3 taken 48294905 times.
✗ Branch 4 not taken.
48294905 std::vector<felInt> locdofid(iLocDof.size());
121
1/2
✓ Branch 3 taken 48294905 times.
✗ Branch 4 not taken.
48294905 numGlobalDof.resize(iComp.size() * iLocDof.size());
122
123 #ifndef NDEBUG
124
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 48294905 times.
48294905 if ( iUnknown == -1) {
125 std::cout << "iVar: " << iVar << std::endl;
126 FEL_ERROR("This variable isn't an unknown in this problem !\n");
127 }
128
129
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 48294905 times.
48294905 if ((unsigned) iVar > m_variable.size() )
130 FEL_ERROR("Not Enough variable!");
131
132
2/2
✓ Branch 1 taken 98723483 times.
✓ Branch 2 taken 48294905 times.
147018388 for(std::size_t i=0; i<iComp.size(); ++i) {
133
1/2
✗ Branch 3 not taken.
✓ Branch 4 taken 98723483 times.
98723483 if ( iComp[i] >= m_variable[iVar].numComponent()) {
134 std::ostringstream oconv;
135 oconv << "Not enough components for this variable! ";
136 oconv << "At least " << m_variable[iVar].numComponent() + 1u << " were expected but "
137 << iComp[i] << " were found.";
138 FEL_ERROR(oconv.str());
139 }
140 }
141 #endif
142
143 48294905 felInt numDofExisting = 0;
144 48294905 felInt numSupportOfDof = 0;
145 // e.g. for Navier Stokes model, we have two unknowns. In m_pSupportDofUnknown, there are first the velocity dof support and then pressure one.
146 // So if iUnknown==0 (velocity), we directly compute the global dof number. Else if iUnknown==1 (pressure), we have to skip all the velocity global dof numbers (numDofExisting).
147 // TODO : rename "numDofExisting"
148
2/2
✓ Branch 0 taken 18599240 times.
✓ Branch 1 taken 48294905 times.
66894145 for ( int i = 0; i < iUnknown; i++) {
149 18599240 iVar2 = m_unknown.idVariable(i);
150 18599240 numSupportOfDof = m_pSupportDofUnknown[i]->numSupportDof();
151 18599240 numDofExisting += numSupportOfDof * m_variable[iVar2].numComponent() ;
152 }
153
154 // dof are arrange in a 1d std::vector like this: iLocDof[0] - iComp[0], iLocDof[1] - iComp[0], ....
155 48294905 std::size_t cnt = 0;
156
2/2
✓ Branch 1 taken 185790699 times.
✓ Branch 2 taken 48294905 times.
234085604 for(std::size_t i=0; i<iLocDof.size(); ++i)
157 185790699 locdofid[i] = numDofExisting + m_pSupportDofUnknown[iUnknown]->iSupportDof()[m_pSupportDofUnknown[iUnknown]->iEle()[iEle] + iLocDof[i] ]*m_variable[iVar].numComponent();
158
159
2/2
✓ Branch 1 taken 98723483 times.
✓ Branch 2 taken 48294905 times.
147018388 for(std::size_t j=0; j<iComp.size(); ++j) {
160
2/2
✓ Branch 1 taken 385655521 times.
✓ Branch 2 taken 98723483 times.
484379004 for(std::size_t i=0; i<iLocDof.size(); ++i) {
161 385655521 numGlobalDof[cnt] = locdofid[i] + iComp[j];
162 385655521 cnt += 1;
163 }
164 }
165 48294905 }
166
167 /*!
168 \brief Function to get global number of dof from one number of support dof.
169 \param[in] idSupportDof number of the support dof.
170 \param[out] idDof std::vector of dofs.
171 */
172 void Dof::supportDofToDof(felInt idSupportDof, std::vector<felInt>& idDof) const {
173 int idVar;
174 int idVariable;
175 felInt numDofExisting = 0;
176 felInt numSupportOfDof = 0;
177 for (unsigned int iUnknown = 0; iUnknown < m_unknown.size(); iUnknown++) {
178 idVar = m_unknown.idVariable(iUnknown);
179 numDofExisting = 0;
180 numSupportOfDof = 0;
181
182 if(FelisceParam::instance().FusionDof) {
183 idSupportDof = m_pSupportDofUnknown[iUnknown]->listPerm()[idSupportDof];
184 }
185
186 for (unsigned int i = 0; i < iUnknown; i++) {
187 numSupportOfDof = m_pSupportDofUnknown[i]->numSupportDof();
188 idVariable = m_unknown.idVariable(i);
189 numDofExisting += numSupportOfDof * m_variable[idVariable].numComponent();
190 }
191
192 for (std::size_t iComp = 0; iComp < m_variable[idVar].numComponent(); iComp++)
193 idDof.push_back(numDofExisting + idSupportDof*m_variable[idVar].numComponent() + iComp);
194 }
195 }
196
197 /*!
198 \brief Find the global number of the component in the mask array.
199 Useful to identify interaction between variable.
200 \param[in] iUnknown Number of the unknown
201 \param[in] iComp Number of the research component for this unknown
202 \return The global number of the component.
203 */
204 677419853 int Dof::getNumGlobComp(std::size_t iUnknown, int iComp) const {
205 677419853 int result = 0;
206 677419853 int idVar = -1;
207
2/4
✗ Branch 1 not taken.
✓ Branch 2 taken 677419853 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 677419853 times.
677419853 FEL_ASSERT( iUnknown < m_unknown.size() && "Not enough unknown in the system.");
208
2/2
✓ Branch 0 taken 107559032 times.
✓ Branch 1 taken 677419853 times.
784978885 for ( std::size_t iUnk = 0; iUnk < iUnknown; iUnk++) {
209 107559032 idVar = m_unknown.idVariable(iUnk);
210 107559032 result += m_variable[idVar].numComponent();
211 }
212 677419853 result +=iComp;
213 677419853 return result;
214 }
215
216 /*!
217 \brief Find the global number of the dof associated to the support dof of a given variable.
218 \param[in] iVar id of the Variable
219 \param[in] listSuppDof sopport of Dofs
220 \param[out] listDofAssociated the global number of the dofs .
221 */
222 void Dof::identifyDofBySupport(int iVar, std::vector<felInt>& listSuppDof, std::vector<felInt>& listDofAssociated) const {
223
224 felInt totNumDof = m_variable[iVar].numComponent()*listSuppDof.size();
225 listDofAssociated.resize(totNumDof,-1);
226 std::vector<felInt> idDofBySupport;
227 felInt indiceStart[m_variable.size()];
228 indiceStart[0] = 0;
229
230 for (int jv = 1; jv < iVar+1 ; jv++) {
231 indiceStart[jv] = indiceStart[jv-1] + m_variable[jv-1].numComponent();
232 }
233
234 felInt NlistSuppDof = static_cast<felInt>(listSuppDof.size());
235 for (felInt iSpdf = 0; iSpdf < NlistSuppDof; iSpdf++) {
236 supportDofToDof(listSuppDof[iSpdf],idDofBySupport);
237 for (std::size_t iComp = 0; iComp < m_variable[iVar].numComponent(); iComp++) {
238 listDofAssociated[iSpdf*m_variable[iVar].numComponent()+iComp] = idDofBySupport[indiceStart[iVar]+iComp];
239 }
240 }
241 }
242
243 // /*!
244 // \brief Method to evaluate the pattern of the matrix.
245 // We realize an "assembly" loop to browse all dof defined (with all unknowns and components)
246 // and we connect dof each others.
247 // After, we store all this information in 2 arrays ( m_iCSR, m_jCSR) as CSR format.
248 // We are using STL unordered_set to sort dof number and eliminate copies.
249 // \todo Optimize loop organisation.
250 // \todo Benoit, 29/10/13 : remove this function ? It seems to not be used anywhere.
251 // */
252 // void Dof::buildPattern() {
253 // felInt numSupportOfDof = 0;
254 // int idVar1 = -1;
255 // int idVar2 = -1;
256
257 // // Eval number total of dof and store total of dof by unknown in STL vector
258 // m_numDofPerUnknown.resize(m_unknown.size());
259 // for (std::size_t iUnknown = 0; iUnknown < m_unknown.size(); iUnknown++) {
260 // idVar1 = m_unknown.idVariable(iUnknown);
261 // numSupportOfDof = m_pSupportDofUnknown[iUnknown]->numSupportDof();
262 // m_numDofPerUnknown[iUnknown] = numSupportOfDof * m_variable[idVar1].numComponent();
263 // m_numDof += m_numDofPerUnknown[iUnknown];
264 // }
265
266 // // loop to connect all dof each others.
267 // std::vector< std::unordered_set<felInt> > nodesNeighborhood(m_numDof);
268 // felInt node1 = 0;
269 // felInt node2 = 0;
270 // for ( std::size_t iUnknown = 0; iUnknown < m_unknown.size(); iUnknown++) {
271 // idVar1 = m_unknown.idVariable(iUnknown);
272 // for (felInt iel = 0; iel < static_cast<felInt>(m_pSupportDofUnknown[iUnknown]->iEle().size()) - 1; iel++) {
273 // for ( int iSupport = 0; iSupport < m_pSupportDofUnknown[iUnknown]->getNumSupportDof(iel); iSupport++) {
274 // for (std::size_t iComp = 0; iComp < m_variable[idVar1].numComponent(); iComp++) {
275 // this->loc2glob(iel,iSupport,idVar1,iComp,node1);
276 // for ( std::size_t iUnknown2 = 0; iUnknown2 < m_unknown.size(); iUnknown2++) {
277 // idVar2 = m_unknown.idVariable(iUnknown2);
278 // for ( int jSupport = 0; jSupport < m_pSupportDofUnknown[iUnknown2]->getNumSupportDof(iel); jSupport++) {
279 // for (std::size_t jComp = 0; jComp < m_variable[idVar2].numComponent(); jComp++) {
280 // const int iConnect = getNumGlobComp( iUnknown, iComp);
281 // const int jConnect = getNumGlobComp( iUnknown2, jComp);
282 // if ( m_unknown.mask()(iConnect, jConnect) > 0) {
283 // this->loc2glob(iel,jSupport,idVar2,jComp,node2);
284 // nodesNeighborhood[node1].insert(node2);
285 // nodesNeighborhood[node2].insert(node1);
286 // }
287 // }
288 // }
289 // }
290 // }
291 // }
292 // }
293 // }
294
295 // //Storage in CSR format
296 // felInt dofSize = 0;
297 // for (unsigned int iNode = 0; iNode < nodesNeighborhood.size(); iNode++) {
298 // dofSize += nodesNeighborhood[iNode].size();
299 // }
300 // m_pattern.resize(m_numDof,dofSize);
301 // felInt pos = 0;
302 // for (unsigned int iNode = 0; iNode < nodesNeighborhood.size(); iNode++) {
303 // m_pattern.rowPointer(iNode+1) = m_pattern.rowPointer(iNode) + nodesNeighborhood[iNode].size();
304 // pos = 0;
305 // for(auto iCon = nodesNeighborhood[iNode].begin(); iCon != nodesNeighborhood[iNode].end(); iCon++) {
306 // m_pattern.columnIndex(m_pattern.rowPointer(iNode) + pos) = *iCon;
307 // pos++;
308 // }
309 // }
310 // }
311
312 /*!
313 \brief Method to evaluate the pattern of the matrix.
314 We realize an "assembly" loop to browse all dof defined (with all unknowns and components)
315 and we connect dof each others.
316 After, we store all this information in 2 arrays ( m_iCSR, m_jCSR) as CSR format.
317 We are using STL unordered_set to sort dof number and eliminate copies.
318 \todo Optimize loop organisation.
319 */
320 525 void Dof::buildPattern(int rankProc, const felInt* dofRepartition) {
321 525 m_pattern.clear();
322
323 525 felInt numDofPerProc = 0;
324
2/2
✓ Branch 0 taken 2608384 times.
✓ Branch 1 taken 525 times.
2608909 for (felInt iDof = 0; iDof < m_numDof; iDof++) {
325
2/2
✓ Branch 0 taken 754597 times.
✓ Branch 1 taken 1853787 times.
2608384 if (dofRepartition[iDof] == rankProc)
326 754597 numDofPerProc++;
327 }
328
329 int idVar1, idVar2;
330 int iMehs1, iMesh2;
331
1/2
✓ Branch 2 taken 525 times.
✗ Branch 3 not taken.
525 std::vector< std::unordered_set<felInt> > nodesNeighborhood(m_numDof);
332 525 felInt node1 = 0;
333 525 felInt node2 = 0;
334
335 // Check
336 525 const auto size_unknown = m_unknown.size();
337
2/4
✗ Branch 1 not taken.
✓ Branch 2 taken 525 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 525 times.
525 FEL_ASSERT(m_pSupportDofUnknown.size() == size_unknown && "Inconsistent size of the system");
338
339 // Loop to connect all dof each others.
340 // For all unknown of the linear system
341
2/2
✓ Branch 0 taken 762 times.
✓ Branch 1 taken 525 times.
1287 for ( std::size_t iUnknown1 = 0; iUnknown1 < size_unknown; iUnknown1++) {
342 762 idVar1 = m_unknown.idVariable(iUnknown1);
343 762 iMehs1 = m_variable[idVar1].idMesh();
344
345 // For all support element of the current unknown
346 762 const auto& r_support_1 = m_pSupportDofUnknown[iUnknown1];
347 762 const auto& r_variable_1 = m_variable[idVar1];
348
2/2
✓ Branch 2 taken 5100262 times.
✓ Branch 3 taken 762 times.
5101024 for (felInt iel = 0; iel < static_cast<felInt>(r_support_1->iEle().size()) - 1; iel++) {
349 // For all support dof of the current support element
350
2/2
✓ Branch 1 taken 19762060 times.
✓ Branch 2 taken 5100262 times.
24862322 for ( int iSupport = 0; iSupport < r_support_1->getNumSupportDof(iel); iSupport++) {
351
352 // For all component of the unknown
353
2/2
✓ Branch 1 taken 42559871 times.
✓ Branch 2 taken 19762060 times.
62321931 for (std::size_t iComp = 0; iComp < r_variable_1.numComponent(); iComp++) {
354
355 // Get the global number of the support dof
356
1/2
✓ Branch 1 taken 42559871 times.
✗ Branch 2 not taken.
42559871 this->loc2glob(iel, iSupport, idVar1, iComp, node1);
357
358 // FEL_ASSERT(m_numDof > node1 && "Inconsistent size of the system");
359
360 // If the support dof is owned by the process
361
2/2
✓ Branch 0 taken 11873805 times.
✓ Branch 1 taken 30686066 times.
42559871 if (dofRepartition[node1] == rankProc) {
362
363 // for all unknown of the linear system
364
2/2
✓ Branch 0 taken 18334438 times.
✓ Branch 1 taken 11873805 times.
30208243 for ( std::size_t iUnknown2 = 0; iUnknown2 < size_unknown; iUnknown2++) {
365 18334438 idVar2 = m_unknown.idVariable(iUnknown2);
366 18334438 iMesh2 = m_variable[idVar2].idMesh();
367
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 18334438 times.
18334438 if ( iMehs1 != iMesh2 )
368 continue;
369
370 // for all support dof of the current support element
371 18334438 const auto& r_support_2 = m_pSupportDofUnknown[iUnknown2];
372 18334438 const auto& r_variable_2 = m_variable[idVar2];
373
2/2
✓ Branch 1 taken 73456036 times.
✓ Branch 2 taken 18334438 times.
91790474 for ( int jSupport = 0; jSupport < r_support_2->getNumSupportDof(iel); jSupport++) {
374
375 // for all component of the second unknown
376
2/2
✓ Branch 1 taken 163816295 times.
✓ Branch 2 taken 73456036 times.
237272331 for (std::size_t jComp = 0; jComp < r_variable_2.numComponent(); jComp++) {
377
378
1/2
✓ Branch 1 taken 163816295 times.
✗ Branch 2 not taken.
163816295 const int iConnect = getNumGlobComp( iUnknown1, iComp);
379
1/2
✓ Branch 1 taken 163816295 times.
✗ Branch 2 not taken.
163816295 const int jConnect = getNumGlobComp( iUnknown2, jComp);
380
381 // FEL_ASSERT(int(m_unknown.mask().size1()) > iConnect && "Inconsistent size of the system");
382 // FEL_ASSERT(int(m_unknown.mask().size2()) > jConnect && "Inconsistent size of the system");
383
384 // if the two components are connected to each other
385
2/4
✓ Branch 2 taken 163816295 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 163816295 times.
✗ Branch 5 not taken.
163816295 if ( m_unknown.mask()(iConnect, jConnect) > 0) {
386 // get the global number of the second support dof
387
1/2
✓ Branch 1 taken 163816295 times.
✗ Branch 2 not taken.
163816295 this->loc2glob(iel,jSupport,idVar2,jComp,node2);
388
389 // FEL_ASSERT(m_numDof > node2 && "Inconsistent size of the system");
390
391 // and insert it in nodesNeighboorhood if it's not in yet.
392
1/2
✓ Branch 2 taken 163816295 times.
✗ Branch 3 not taken.
163816295 nodesNeighborhood[node1].insert(node2); // TODO D.C. why we don't check if ( node1 != node2) as everywhere else???
393 }
394 }
395 }
396 }
397 }
398 }
399 }
400 }
401 }
402
403 //Storage in CSR format
404 525 felInt dofSize = 0;
405
2/2
✓ Branch 1 taken 2608384 times.
✓ Branch 2 taken 525 times.
2608909 for (unsigned int iNode = 0; iNode < nodesNeighborhood.size(); iNode++) {
406 2608384 dofSize += nodesNeighborhood[iNode].size();
407 }
408
1/2
✓ Branch 1 taken 525 times.
✗ Branch 2 not taken.
525 m_pattern.resize(numDofPerProc,dofSize);
409 525 felInt pos = 0;
410 525 felInt cptDof = 0;
411
2/2
✓ Branch 1 taken 2608384 times.
✓ Branch 2 taken 525 times.
2608909 for (unsigned int iNode = 0; iNode < nodesNeighborhood.size(); iNode++) {
412
2/2
✓ Branch 0 taken 754597 times.
✓ Branch 1 taken 1853787 times.
2608384 if (dofRepartition[iNode] == rankProc) {
413
2/4
✓ Branch 1 taken 754597 times.
✗ Branch 2 not taken.
✓ Branch 6 taken 754597 times.
✗ Branch 7 not taken.
754597 m_pattern.rowPointer(cptDof+1) = m_pattern.rowPointer(cptDof) + nodesNeighborhood[iNode].size();
414 754597 pos = 0;
415
2/2
✓ Branch 5 taken 31845273 times.
✓ Branch 6 taken 754597 times.
32599870 for(auto iCon = nodesNeighborhood[iNode].begin(); iCon != nodesNeighborhood[iNode].end(); iCon++) {
416
2/4
✓ Branch 2 taken 31845273 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 31845273 times.
✗ Branch 6 not taken.
31845273 m_pattern.columnIndex(m_pattern.rowPointer(cptDof) + pos) = *iCon;
417 31845273 pos++;
418 }
419 754597 cptDof++;
420 }
421 }
422 525 }
423
424 /*!
425 \brief Method to initialize the pattern with a uniform repartition of the dof on the processes.
426 We realize an "assembly" loop to browse all dof defined (with all unknowns and components)
427 and we connect dof each others.
428 After, we store all this information in 2 arrays ( m_iCSR, m_jCSR) as CSR format.
429 We are using STL set to sort dof number and eliminate copies.
430 \todo Optimize loop organisation.
431 */
432 525 void Dof::initializePattern(int sizeProc, int rankProc) {
433 525 felInt numSupportOfDof = 0;
434 int idVar1, idVar2;
435 int iMesh1, iMesh2;
436 // Eval number total of dof and store total of dof by unknown in STL vector
437 525 const auto size_unknown = m_unknown.size();
438
1/2
✓ Branch 1 taken 525 times.
✗ Branch 2 not taken.
525 m_numDofPerUnknown.resize(size_unknown);
439
440
2/2
✓ Branch 0 taken 762 times.
✓ Branch 1 taken 525 times.
1287 for (std::size_t iUnknown = 0; iUnknown < size_unknown; iUnknown++) {
441 762 idVar1 = m_unknown.idVariable(iUnknown);
442 762 numSupportOfDof = m_pSupportDofUnknown[iUnknown]->numSupportDof();
443 762 m_numDofPerUnknown[iUnknown] = numSupportOfDof * m_variable[idVar1].numComponent();
444 762 m_numDof += m_numDofPerUnknown[iUnknown];
445 }
446
447 // number of dof in this proc
448 525 felInt numDofByProc = m_numDof/sizeProc;
449 525 felInt beginIdDofLocal = rankProc*numDofByProc;
450 felInt numDofProc;
451 felInt endIdDofLocal;
452
2/2
✓ Branch 0 taken 155 times.
✓ Branch 1 taken 370 times.
525 if(rankProc == sizeProc-1) {
453 155 numDofProc = m_numDof - beginIdDofLocal;
454 155 endIdDofLocal = m_numDof;
455 } else {
456 370 numDofProc = numDofByProc;
457 370 endIdDofLocal = beginIdDofLocal + numDofProc;
458 }
459
460 // Loop to connect all dof each others.
461
1/2
✓ Branch 2 taken 525 times.
✗ Branch 3 not taken.
525 std::vector< std::unordered_set<felInt> > nodesNeighborhood(numDofProc);
462 525 felInt node1 = 0;
463 525 felInt node2 = 0;
464
465 // For all unknowns
466
2/2
✓ Branch 0 taken 762 times.
✓ Branch 1 taken 525 times.
1287 for ( std::size_t iUnknown1 = 0; iUnknown1 < size_unknown; iUnknown1++) {
467 762 idVar1 = m_unknown.idVariable(iUnknown1);
468 762 iMesh1 = m_variable[idVar1].idMesh();
469
470 // For all global support element of this unknown
471 762 const auto& r_support_1 = m_pSupportDofUnknown[iUnknown1];
472 762 const auto& r_variable_1 = m_variable[idVar1];
473
2/2
✓ Branch 2 taken 5100262 times.
✓ Branch 3 taken 762 times.
5101024 for (felInt iel = 0; iel < static_cast<felInt>(r_support_1->iEle().size()) - 1; iel++) {
474
475 // For all support dof in this element
476
2/2
✓ Branch 1 taken 19762060 times.
✓ Branch 2 taken 5100262 times.
24862322 for ( int iSupport = 0; iSupport < r_support_1->getNumSupportDof(iel); iSupport++) {
477
478 // For all component of the unknown
479
2/2
✓ Branch 1 taken 42559871 times.
✓ Branch 2 taken 19762060 times.
62321931 for (std::size_t iComp = 0; iComp < r_variable_1.numComponent(); iComp++) {
480 // Get the global id of the support dof
481
1/2
✓ Branch 1 taken 42559871 times.
✗ Branch 2 not taken.
42559871 loc2glob(iel, iSupport, idVar1, iComp, node1);
482
483 // If this support dof is owned by this process
484
4/4
✓ Branch 0 taken 27391755 times.
✓ Branch 1 taken 15168116 times.
✓ Branch 2 taken 11873805 times.
✓ Branch 3 taken 15517950 times.
42559871 if(node1 >= beginIdDofLocal && node1 < endIdDofLocal) {
485 // For all unknown in the linear problem
486
2/2
✓ Branch 0 taken 18334438 times.
✓ Branch 1 taken 11873805 times.
30208243 for ( std::size_t iUnknown2 = 0; iUnknown2 < size_unknown; iUnknown2++) {
487 18334438 idVar2 = m_unknown.idVariable(iUnknown2);
488 18334438 iMesh2 = m_variable[idVar2].idMesh();
489
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 18334438 times.
18334438 if ( iMesh1 != iMesh2 )
490 continue;
491
492 // For all support dof in the current element
493 18334438 const auto& r_support_2 = m_pSupportDofUnknown[iUnknown2];
494 18334438 const auto& r_variable_2 = m_variable[idVar2];
495
2/2
✓ Branch 1 taken 73456036 times.
✓ Branch 2 taken 18334438 times.
91790474 for (int jSupport = 0; jSupport < r_support_2->getNumSupportDof(iel); jSupport++) {
496
497 // For all component of the second unknown
498
2/2
✓ Branch 1 taken 163816295 times.
✓ Branch 2 taken 73456036 times.
237272331 for (std::size_t jComp = 0; jComp < r_variable_2.numComponent(); jComp++) {
499 // If the two current components are connected
500
1/2
✓ Branch 1 taken 163816295 times.
✗ Branch 2 not taken.
163816295 const int iConnect = getNumGlobComp( iUnknown1, iComp);
501
1/2
✓ Branch 1 taken 163816295 times.
✗ Branch 2 not taken.
163816295 const int jConnect = getNumGlobComp( iUnknown2, jComp);
502
2/4
✓ Branch 2 taken 163816295 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 163816295 times.
✗ Branch 5 not taken.
163816295 if ( m_unknown.mask()(iConnect, jConnect) > 0) {
503 // get the global id of the second support dof
504
1/2
✓ Branch 1 taken 163816295 times.
✗ Branch 2 not taken.
163816295 loc2glob(iel, jSupport, idVar2, jComp, node2);
505
506 // remove diagonal term to define pattern and use it in Parmetis.
507 // if the two global ids are different, they are neighboors
508
2/2
✓ Branch 0 taken 151942238 times.
✓ Branch 1 taken 11874057 times.
163816295 if ( node1 != node2)
509
1/2
✓ Branch 2 taken 151942238 times.
✗ Branch 3 not taken.
151942238 nodesNeighborhood[node1-beginIdDofLocal].insert(node2);
510 }
511 }
512 }
513 }
514 }
515 }
516 }
517 }
518 }
519
520 // Storage in CSR format
521 525 std::size_t dofSize = 0;
522
2/2
✓ Branch 1 taken 754597 times.
✓ Branch 2 taken 525 times.
755122 for (std::size_t iNode = 0; iNode < nodesNeighborhood.size(); iNode++) {
523 754597 dofSize += nodesNeighborhood[iNode].size();
524 }
525
1/2
✓ Branch 1 taken 525 times.
✗ Branch 2 not taken.
525 m_pattern.resize(numDofProc,dofSize);
526 525 felInt pos = 0;
527
2/2
✓ Branch 0 taken 754597 times.
✓ Branch 1 taken 525 times.
755122 for (felInt iNode = 0; iNode < numDofProc; iNode++) {
528
2/4
✓ Branch 1 taken 754597 times.
✗ Branch 2 not taken.
✓ Branch 6 taken 754597 times.
✗ Branch 7 not taken.
754597 m_pattern.rowPointer(iNode+1) = m_pattern.rowPointer(iNode) + nodesNeighborhood[iNode].size();
529 754597 pos = 0;
530
2/2
✓ Branch 5 taken 31090676 times.
✓ Branch 6 taken 754597 times.
31845273 for(auto iCon = nodesNeighborhood[iNode].begin(); iCon != nodesNeighborhood[iNode].end(); iCon++) {
531
2/4
✓ Branch 2 taken 31090676 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 31090676 times.
✗ Branch 6 not taken.
31090676 m_pattern.columnIndex(m_pattern.rowPointer(iNode) + pos) = *iCon;
532 31090676 pos++;
533 }
534 }
535 525 }
536
537 21 void Dof::mergeGlobalPattern(const std::vector<felInt>& iCSRc_new, const std::vector<felInt>& jCSRc_new) {
538 // Size of iCSR (number of dof + 1)
539 21 std::size_t dimICSR = m_pattern.numRows() + 1;
540
541 // Check that the given complementary pattern has the right size
542
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 21 times.
21 FEL_ASSERT(iCSRc_new.size() == dimICSR);
543
544 // if m_patternC is not created already, initialize it to nothing
545
1/2
✓ Branch 0 taken 21 times.
✗ Branch 1 not taken.
21 if(m_patternIsChanged == 0)
546
1/2
✓ Branch 1 taken 21 times.
✗ Branch 2 not taken.
21 m_patternC.resize(dimICSR, 0);
547
548 // temporary set to add the new non zero element contains in the new complementary pattern
549 21 std::set<felInt> tmpset, tmpsetc;
550 21 std::pair<std::set<felInt>::iterator, bool> ret;
551
552 // New jCSR and iCSR for m_pattern and m_patternC
553
2/4
✓ Branch 2 taken 21 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 21 times.
✗ Branch 7 not taken.
42 std::vector<felInt> iCSR(dimICSR), iCSRc(dimICSR);
554 21 std::vector<felInt> jCSR, jCSRc;
555
556 // Initialize first entry of iCSR and iCSRc
557 21 iCSR[0] = 0;
558 21 iCSRc[0] = 0;
559
560 // Reserve space for jCSR and jCSRc
561
1/2
✓ Branch 3 taken 21 times.
✗ Branch 4 not taken.
21 jCSR.reserve(m_pattern.numNonzeros() + jCSRc_new.size());
562
1/2
✓ Branch 3 taken 21 times.
✗ Branch 4 not taken.
21 jCSRc.reserve(m_patternC.numNonzeros() + jCSRc_new.size());
563
564 // loop over each dof
565
2/2
✓ Branch 0 taken 73181 times.
✓ Branch 1 taken 21 times.
73202 for(std::size_t n=0; n<dimICSR-1; ++n) {
566 // fill tmpset with the current indices of column for this row
567 73181 tmpset.clear();
568
4/6
✓ Branch 1 taken 73181 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1493501 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1420320 times.
✓ Branch 7 taken 73181 times.
1493501 for(felInt i=m_pattern.rowPointer(n); i<m_pattern.rowPointer(n+1); ++i)
569
2/4
✓ Branch 1 taken 1420320 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1420320 times.
✗ Branch 5 not taken.
1420320 tmpset.insert(m_pattern.columnIndex(i));
570
571 // fill tmpsetc with the current indices of column for this row
572 73181 tmpsetc.clear();
573
3/6
✓ Branch 1 taken 73181 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 73181 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 73181 times.
73181 for(felInt i=m_patternC.rowPointer(n); i<m_patternC.rowPointer(n+1); ++i)
574 tmpsetc.insert(m_patternC.columnIndex(i));
575
576 // insert the new nonzero of jCSRc_new in the sets
577
2/2
✓ Branch 2 taken 2666302 times.
✓ Branch 3 taken 73181 times.
2739483 for(felInt i=iCSRc_new[n]; i<iCSRc_new[n+1]; ++i) {
578
1/2
✓ Branch 2 taken 2666302 times.
✗ Branch 3 not taken.
2666302 ret = tmpset.insert(jCSRc_new[i]);
579
2/2
✓ Branch 0 taken 1246138 times.
✓ Branch 1 taken 1420164 times.
2666302 if(ret.second)
580
1/2
✓ Branch 2 taken 1246138 times.
✗ Branch 3 not taken.
1246138 tmpsetc.insert(jCSRc_new[i]);
581 }
582
583 // set the iCSR component of the patterns
584 73181 iCSR[n+1] = iCSR[n] + tmpset.size();
585 73181 iCSRc[n+1] = iCSRc[n] + tmpsetc.size();
586
587 // set the jCSR component of the patterns
588
2/2
✓ Branch 4 taken 2666458 times.
✓ Branch 5 taken 73181 times.
2739639 for (auto it = tmpset.begin(); it != tmpset.end(); ++it)
589
1/2
✓ Branch 2 taken 2666458 times.
✗ Branch 3 not taken.
2666458 jCSR.push_back(*it);
590
591
2/2
✓ Branch 4 taken 1246138 times.
✓ Branch 5 taken 73181 times.
1319319 for (auto it = tmpsetc.begin(); it != tmpsetc.end(); ++it)
592
1/2
✓ Branch 2 taken 1246138 times.
✗ Branch 3 not taken.
1246138 jCSRc.push_back(*it);
593 }
594
595 // Resize the jCSR componenent of the patterns (remove the not needed reserve space)
596
1/2
✓ Branch 2 taken 21 times.
✗ Branch 3 not taken.
21 jCSR.resize(iCSR[dimICSR-1]);
597
1/2
✓ Branch 2 taken 21 times.
✗ Branch 3 not taken.
21 jCSRc.resize(iCSRc[dimICSR-1]);
598
599 // Set the new patterns
600
1/2
✓ Branch 1 taken 21 times.
✗ Branch 2 not taken.
21 m_pattern.set(iCSR, jCSR);
601
1/2
✓ Branch 1 taken 21 times.
✗ Branch 2 not taken.
21 m_patternC.set(iCSRc, jCSRc);
602
603 // Set the flag indicating that there is a non empty complementary pattern
604 21 m_patternIsChanged = true;
605 21 }
606
607 /*!
608 \brief Build complementary pattern with the new repartition of the dof
609 \param[in] mpiComm MPI communicator
610 \param[in] numProc Total number of process.
611 \param[in] rankProc Rank of the current process.
612 \param[in] dofRepartition New repartition of the dof among the process.
613 */
614 21 void Dof::mergeLocalCompPattern(MPI_Comm mpiComm, felInt numProc, felInt rankProc, const std::vector<felInt>& dofRepartition) {
615 // build the uniform initial repartition
616 21 felInt numDofByProc = m_numDof/numProc;
617 21 felInt beginIdDofLocal = rankProc*numDofByProc;
618 felInt endIdDofLocal;
619 felInt numDofProc;
620
2/2
✓ Branch 0 taken 15 times.
✓ Branch 1 taken 6 times.
21 if(rankProc == numProc-1) {
621 15 numDofProc = m_numDof - beginIdDofLocal;
622 15 endIdDofLocal = m_numDof;
623 } else {
624 6 numDofProc = numDofByProc;
625 6 endIdDofLocal = beginIdDofLocal + numDofProc;
626 }
627
628 // Complementary pattern on the new repartition of dof
629
1/2
✓ Branch 3 taken 21 times.
✗ Branch 4 not taken.
42 std::vector<felInt> iCSRc(m_pattern.numRows() + 1), jCSRc;
630
631 // New pattern
632 21 std::vector<felInt> iCSR, jCSR;
633
634 // Initialize first value of iCSRc and reserve space for jCSRc
635 21 iCSRc[0] = 0;
636
1/2
✓ Branch 2 taken 21 times.
✗ Branch 3 not taken.
21 jCSRc.reserve(m_patternC.numNonzeros());
637
638 // Count the dof own by this process
639 21 felInt countDof = 0;
640 21 felInt countDofInit = 0;
641 21 felInt numNonZeros = 0;
642
643 // Temporary vector
644 21 std::vector<felInt> tmpvec;
645
646 // MPI status
647 MPI_Status mpiStatus;
648
649 // build new complementary pattern
650
2/2
✓ Branch 0 taken 76589 times.
✓ Branch 1 taken 21 times.
76610 for(felInt idof=0; idof<m_numDof; ++idof) {
651 // if the dof is own by this process
652
2/2
✓ Branch 1 taken 73181 times.
✓ Branch 2 taken 3408 times.
76589 if(dofRepartition[idof] == rankProc) {
653 // if the dof was also on this process with the initial repartition
654
4/4
✓ Branch 0 taken 72762 times.
✓ Branch 1 taken 419 times.
✓ Branch 2 taken 72333 times.
✓ Branch 3 taken 429 times.
73181 if((idof >= beginIdDofLocal) && (idof < endIdDofLocal)) {
655 // take directly from the current complementary pattern
656
1/2
✓ Branch 1 taken 72333 times.
✗ Branch 2 not taken.
72333 numNonZeros = m_patternC.numNonzerosInRow(countDofInit);
657 72333 iCSRc[countDof + 1] = iCSRc[countDof] + numNonZeros;
658
4/6
✓ Branch 1 taken 72333 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1298509 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1226176 times.
✓ Branch 7 taken 72333 times.
1298509 for(felInt i=m_patternC.rowPointer(countDofInit); i<m_patternC.rowPointer(countDofInit+1); ++i)
659
2/4
✓ Branch 1 taken 1226176 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1226176 times.
✗ Branch 5 not taken.
1226176 jCSRc.push_back(m_patternC.columnIndex(i));
660
661 72333 ++countDofInit;
662 72333 ++countDof;
663 72333 } else {
664 // ask the process owning it in the initial repartition
665 848 felInt senderRank = 1;
666
2/2
✓ Branch 0 taken 1277 times.
✓ Branch 1 taken 848 times.
2125 while(idof >= senderRank * numDofByProc)
667 1277 ++senderRank;
668 848 senderRank = std::min(senderRank - 1, numProc - 1);
669
670
1/2
✓ Branch 1 taken 848 times.
✗ Branch 2 not taken.
848 MPI_Probe(senderRank, idof, mpiComm, &mpiStatus);
671
1/2
✓ Branch 1 taken 848 times.
✗ Branch 2 not taken.
848 MPI_Get_count(&mpiStatus, MPI_INT, &numNonZeros);
672
1/2
✓ Branch 1 taken 848 times.
✗ Branch 2 not taken.
848 tmpvec.resize(numNonZeros);
673
1/2
✓ Branch 2 taken 848 times.
✗ Branch 3 not taken.
848 MPI_Recv(tmpvec.data(), numNonZeros, MPI_INT, senderRank, idof, mpiComm, &mpiStatus);
674
2/2
✓ Branch 1 taken 19962 times.
✓ Branch 2 taken 848 times.
20810 for(std::size_t i=0; i<tmpvec.size(); ++i)
675
1/2
✓ Branch 2 taken 19962 times.
✗ Branch 3 not taken.
19962 jCSRc.push_back(tmpvec[i]);
676
677 848 iCSRc[countDof + 1] = iCSRc[countDof] + numNonZeros;
678 848 ++countDof;
679 }
680 } else {
681 // if this dof is own by this process with the initial repartition
682
4/4
✓ Branch 0 taken 2123 times.
✓ Branch 1 taken 1285 times.
✓ Branch 2 taken 848 times.
✓ Branch 3 taken 1275 times.
3408 if((idof >= beginIdDofLocal) && (idof < endIdDofLocal)) {
683 // send to the new process
684 848 felInt localPos = idof - beginIdDofLocal;
685
1/2
✓ Branch 1 taken 848 times.
✗ Branch 2 not taken.
848 numNonZeros = m_patternC.numNonzerosInRow(localPos);
686
687 // We, Benoit and Matteo, added this check. (27 may 2015)
688 // The reason behind is that, even when sending zero integers,
689 // the address of "m_patternC.columnIndex(m_patternC.rowPointer(localPos))" was
690 // evaluated and led to a segmentation fault when accessing to m_colind into the complementary pattern
691 // ( at least in debug mode, in release it was not clear)
692 // This happened, for instance, when the complementary pattern was completely empty on one processor (think about implicit windkessel for instance).
693 // It was not possible to avoid at all the send when sending zero integers since the other processor would have waited in any case:
694
1/2
✓ Branch 0 taken 848 times.
✗ Branch 1 not taken.
848 if ( numNonZeros > 0 ) {
695
3/6
✓ Branch 2 taken 848 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 848 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 848 times.
✗ Branch 9 not taken.
848 MPI_Send(&m_patternC.columnIndex(m_patternC.rowPointer(localPos)), numNonZeros, MPI_INT, dofRepartition[idof], idof, mpiComm);
696 }
697 else {
698 int trash(-1);
699 MPI_Send(&trash, numNonZeros, MPI_INT, dofRepartition[idof], idof, mpiComm);
700 }
701 848 ++countDofInit;
702 }
703 }
704 }
705
706 // Add the complementary pattern to the pattern
707
1/2
✓ Branch 3 taken 21 times.
✗ Branch 4 not taken.
21 jCSRc.resize(iCSRc[iCSRc.size() - 1]);
708
1/2
✓ Branch 2 taken 21 times.
✗ Branch 3 not taken.
21 iCSR.resize(iCSRc.size());
709
1/2
✓ Branch 3 taken 21 times.
✗ Branch 4 not taken.
21 jCSR.resize(m_pattern.numNonzeros() + jCSRc.size());
710
711 21 iCSR[0] = 0;
712 21 felInt count = 0;
713
2/2
✓ Branch 1 taken 73181 times.
✓ Branch 2 taken 21 times.
73202 for(std::size_t i=0; i<iCSR.size()-1; ++i) {
714 // iCSR
715
1/2
✓ Branch 2 taken 73181 times.
✗ Branch 3 not taken.
73181 iCSR[i+1] = iCSR[i] + m_pattern.numNonzerosInRow(i) + (iCSRc[i+1] - iCSRc[i]);
716
717 // jCSR
718
4/6
✓ Branch 1 taken 73181 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1566682 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1493501 times.
✓ Branch 7 taken 73181 times.
1566682 for(felInt j=m_pattern.rowPointer(i); j<m_pattern.rowPointer(i+1); ++j) {
719
1/2
✓ Branch 1 taken 1493501 times.
✗ Branch 2 not taken.
1493501 jCSR[count] = m_pattern.columnIndex(j);
720 1493501 ++count;
721 }
722
723
2/2
✓ Branch 2 taken 1246138 times.
✓ Branch 3 taken 73181 times.
1319319 for(felInt j=iCSRc[i]; j<iCSRc[i+1]; ++j) {
724 1246138 jCSR[count] = jCSRc[j];
725 1246138 ++count;
726 }
727 }
728
729 // set the patterns
730
1/2
✓ Branch 1 taken 21 times.
✗ Branch 2 not taken.
21 m_patternC.set(iCSRc, jCSRc);
731
1/2
✓ Branch 1 taken 21 times.
✗ Branch 2 not taken.
21 m_pattern.set(iCSR, jCSR);
732 21 }
733
734 void Dof::getNeighbourhood( felInt node, std::vector<felInt> & neighbourhood ) const{
735 neighbourhood.clear();
736 const std::size_t start = m_pattern.rowPointer( node );
737 const std::size_t size = m_pattern.rowPointer( node + 1 ) - start;
738 for ( std::size_t k(0); k < size; ++k)
739 neighbourhood.push_back( m_pattern.columnIndex( start + k ) );
740 }
741
742 /*!
743 * Print pattern of the matrix and number of dofs associate at each variable.
744 */
745 260 void Dof::print(int verbose, std::ostream& outstr) const {
746
1/2
✓ Branch 0 taken 260 times.
✗ Branch 1 not taken.
260 if (verbose) {
747 260 outstr << "Dof: " << m_numDof << std::endl;
748 260 outstr << "Number of dof per unknown:\n";
749
2/2
✓ Branch 1 taken 429 times.
✓ Branch 2 taken 260 times.
689 for (unsigned int iNumDof = 0; iNumDof < m_numDofPerUnknown.size(); iNumDof++)
750 429 outstr << "unknown " << iNumDof << ": " << m_numDofPerUnknown[iNumDof] << "\n";
751 260 outstr << std::endl;
752
753
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 260 times.
260 if(verbose > 25) {
754 outstr << "PATTERN" << std::endl;
755 m_pattern.print(verbose, outstr);
756
757 if(m_patternIsChanged == 1) {
758 outstr << "COMPLEMENTARY PATTERN" << std::endl;
759 m_patternC.print(verbose, outstr);
760 }
761 }
762 }
763 260 }
764 }
765
766