GCC Code Coverage Report


Directory: ./
File: FiniteElement/curBaseFiniteElement.cpp
Date: 2024-04-14 07:32:34
Exec Total Coverage
Lines: 125 193 64.8%
Branches: 94 189 49.7%

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: J-F. Gerbeau and V. Martin
13 //
14
15 // System includes
16 #include <cmath>
17
18 // External includes
19
20 // Project includes
21 #include "FiniteElement/curBaseFiniteElement.hpp"
22 #include "Tools/math_utilities.hpp"
23
24 namespace felisce
25 {
26 170305 CurBaseFiniteElement::CurBaseFiniteElement(
27 const RefElement& refEle,
28 const GeoElement& geoEle,
29 const std::vector<DegreeOfExactness>& degOfExactness
30 170305 ):
31 170305 m_refEle(refEle),
32 170305 m_geoEle(geoEle),
33 170305 m_numDof(m_refEle.numDOF()),
34 170305 m_numRefCoor(m_geoEle.numCoor()),
35 170305 m_numCoor(m_numRefCoor), // will be changed by curvilinearFE
36
3/6
✓ Branch 2 taken 170305 times.
✗ Branch 3 not taken.
✓ Branch 10 taken 170305 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 170305 times.
✗ Branch 14 not taken.
340610 m_numPoint(m_geoEle.numPoint())
37 {
38
7/8
✓ Branch 2 taken 228 times.
✓ Branch 3 taken 87517 times.
✓ Branch 4 taken 56958 times.
✓ Branch 5 taken 8350 times.
✓ Branch 6 taken 9456 times.
✓ Branch 7 taken 39 times.
✓ Branch 8 taken 7757 times.
✗ Branch 9 not taken.
170305 switch (m_geoEle.shape().typeShape())
39 {
40 228 case Node:
41
1/2
✓ Branch 2 taken 228 times.
✗ Branch 3 not taken.
228 m_quadratureRule = &listQuadratureRuleNode.quadratureRuleByExactness(degOfExactness[0]);
42 228 break;
43 87517 case Segment:
44
1/2
✓ Branch 2 taken 87517 times.
✗ Branch 3 not taken.
87517 m_quadratureRule = &listQuadratureRuleSegment.quadratureRuleByExactness(degOfExactness[0]);
45 87517 break;
46 56958 case Triangle:
47 // if (FelisceParam::instance().typeOfShellModel=="MITC3" || FelisceParam::instance().typeOfShellModel=="MITC3+")
48 // m_quadratureRule = &quadratureRuleMitc3;
49 // else
50
1/2
✓ Branch 2 taken 56958 times.
✗ Branch 3 not taken.
56958 m_quadratureRule = &listQuadratureRuleTriangle.quadratureRuleByExactness(degOfExactness[0]);
51 56958 break;
52 8350 case Quadrilateral:
53
1/2
✓ Branch 2 taken 8350 times.
✗ Branch 3 not taken.
8350 m_quadratureRule = &listQuadratureRuleQuadrilateral.quadratureRuleByExactness(degOfExactness[0]);
54 8350 break;
55 9456 case Tetrahedron:
56
1/2
✓ Branch 2 taken 9456 times.
✗ Branch 3 not taken.
9456 m_quadratureRule = &listQuadratureRuleTetrahedron.quadratureRuleByExactness(degOfExactness[0]);
57 9456 break;
58 39 case Hexahedron:
59
2/2
✓ Branch 1 taken 14 times.
✓ Branch 2 taken 25 times.
39 if (degOfExactness.size() == 1) {
60
1/2
✓ Branch 2 taken 14 times.
✗ Branch 3 not taken.
14 m_quadratureRule = &listQuadratureRuleHexahedron.quadratureRuleByExactness(degOfExactness[0]);
61 } else {
62
1/2
✓ Branch 3 taken 25 times.
✗ Branch 4 not taken.
25 m_quadratureRule = &listQuadratureRuleHexahedronCombined.quadratureRuleByExactness(degOfExactness[1], degOfExactness[2]);
63 }
64 39 break;
65 7757 case Prism:
66
2/2
✓ Branch 1 taken 2712 times.
✓ Branch 2 taken 5045 times.
7757 if (degOfExactness.size() == 1) {
67
1/2
✓ Branch 2 taken 2712 times.
✗ Branch 3 not taken.
2712 m_quadratureRule = &listQuadratureRulePrism.quadratureRuleByExactness(degOfExactness[0]);
68 } else {
69
1/2
✓ Branch 3 taken 5045 times.
✗ Branch 4 not taken.
5045 m_quadratureRule = &listQuadratureRulePrismCombined.quadratureRuleByExactness(degOfExactness[1], degOfExactness[2]);
70 }
71 7757 break;
72 default:
73 FEL_ERROR("There is no quadrature rule list defined for this element shape");
74 break;
75 }
76
1/2
✓ Branch 1 taken 170305 times.
✗ Branch 2 not taken.
170305 m_quadratureRule->print(0);
77 // TODO: Replace variable with method and call
78 170305 m_numQuadraturePoint = m_quadratureRule->numQuadraturePoint();
79
80
1/2
✓ Branch 1 taken 170305 times.
✗ Branch 2 not taken.
170305 m_point.resize(m_numPoint,m_numCoor); // can be MODIFIED if used by curvilinearFE
81
1/2
✓ Branch 1 taken 170305 times.
✗ Branch 2 not taken.
170305 currentQuadPoint.resize(m_numQuadraturePoint);
82
1/2
✓ Branch 1 taken 170305 times.
✗ Branch 2 not taken.
170305 m_dPhiRef.resize(m_numQuadraturePoint );
83
1/2
✓ Branch 1 taken 170305 times.
✗ Branch 2 not taken.
170305 m_phiGeo.resize(m_numQuadraturePoint);
84
1/2
✓ Branch 1 taken 170305 times.
✗ Branch 2 not taken.
170305 m_dPhiGeo.resize(m_numQuadraturePoint);
85
1/2
✓ Branch 1 taken 170305 times.
✗ Branch 2 not taken.
170305 m_meas.resize(m_numQuadraturePoint);
86
1/2
✓ Branch 1 taken 170305 times.
✗ Branch 2 not taken.
170305 phi.resize(m_numQuadraturePoint);
87 //dPhi.resize(m_numQuadraturePoint);
88
1/2
✓ Branch 1 taken 170305 times.
✗ Branch 2 not taken.
170305 weightMeas.resize(m_numQuadraturePoint);
89
90
2/2
✓ Branch 0 taken 538635 times.
✓ Branch 1 taken 170305 times.
708940 for(int ig=0; ig<m_numQuadraturePoint; ig++) {
91
1/2
✓ Branch 2 taken 538635 times.
✗ Branch 3 not taken.
538635 m_dPhiRef[ig].resize(m_numRefCoor,m_numDof);
92
1/2
✓ Branch 2 taken 538635 times.
✗ Branch 3 not taken.
538635 m_phiGeo[ig].resize(m_numPoint);
93
1/2
✓ Branch 2 taken 538635 times.
✗ Branch 3 not taken.
538635 m_dPhiGeo[ig].resize(m_numRefCoor,m_numPoint);
94
1/2
✓ Branch 2 taken 538635 times.
✗ Branch 3 not taken.
538635 phi[ig].resize(m_numDof);
95 //dPhi[ig].resize(m_numRefCoor,m_numDof);
96
97
1/2
✓ Branch 1 taken 538635 times.
✗ Branch 2 not taken.
538635 const QuadraturePoint& quadPt = m_quadratureRule->quadraturePoint(ig);
98
2/2
✓ Branch 0 taken 1818696 times.
✓ Branch 1 taken 538635 times.
2357331 for(int i=0; i<m_numDof; i++ ) {
99
2/4
✓ Branch 2 taken 1818696 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 1818696 times.
✗ Branch 7 not taken.
1818696 phi[ig](i) = m_refEle.basisFunction().phi(i, quadPt);
100
2/2
✓ Branch 0 taken 3852792 times.
✓ Branch 1 taken 1818696 times.
5671488 for(int icoor=0; icoor<m_numRefCoor; icoor++ ) {
101
2/4
✓ Branch 2 taken 3852792 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 3852792 times.
✗ Branch 7 not taken.
3852792 m_dPhiRef[ig](icoor,i) = m_refEle.basisFunction().dPhi(i,icoor, quadPt);
102 }
103 }
104
2/2
✓ Branch 0 taken 1755134 times.
✓ Branch 1 taken 538635 times.
2293769 for(int k=0; k<m_numPoint; k++) {
105
2/4
✓ Branch 2 taken 1755134 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 1755134 times.
✗ Branch 7 not taken.
1755134 m_phiGeo[ig](k)= m_geoEle.basisFunction().phi(k, quadPt);
106
2/2
✓ Branch 0 taken 3693076 times.
✓ Branch 1 taken 1755134 times.
5448210 for(int icoor=0; icoor<m_numRefCoor; icoor++) {
107
2/4
✓ Branch 2 taken 3693076 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 3693076 times.
✗ Branch 7 not taken.
3693076 m_dPhiGeo[ig](icoor,k) = m_geoEle.basisFunction().dPhi(k,icoor, quadPt);
108 }
109 }
110 }
111 170305 }
112
113 /***********************************************************************************/
114 /***********************************************************************************/
115
116 117157 CurBaseFiniteElement::CurBaseFiniteElement(
117 const RefElement& refEle,
118 const GeoElement& geoEle,
119 const DegreeOfExactness& degOfExactness
120
2/4
✓ Branch 2 taken 117157 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 117157 times.
✗ Branch 6 not taken.
117157 ) : CurBaseFiniteElement(refEle, geoEle, std::vector<DegreeOfExactness>(1, degOfExactness))
121 {
122 117157 }
123
124 /***********************************************************************************/
125 /***********************************************************************************/
126
127 373212 CurBaseFiniteElement::~CurBaseFiniteElement()
128 {
129 373212 m_point.clear();
130 373212 currentQuadPoint.clear();
131 373212 m_dPhiRef.clear();
132 373212 m_phiGeo.clear();
133 373212 m_dPhiGeo.clear();
134 373212 m_meas.clear();
135 373212 phi.clear();
136 //dPhi.clear();
137 373212 weightMeas.clear();
138 }
139
140 /***********************************************************************************/
141 /***********************************************************************************/
142
143 void CurBaseFiniteElement::placeholder() const
144 {
145 // Do nothing! Exists just to avoid emitting vtables in every translation unit.
146 }
147
148 /***********************************************************************************/
149 /***********************************************************************************/
150
151 764325 void CurBaseFiniteElement::computeCurrentQuadraturePoint()
152 {
153
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 764325 times.
764325 FEL_ASSERT( hasPoint() );
154
155
2/2
✓ Branch 0 taken 3806084 times.
✓ Branch 1 taken 764325 times.
4570409 for (int ig = 0; ig < m_numQuadraturePoint; ig++) {
156 3806084 coorMap(currentQuadPoint[ig],m_quadratureRule->quadraturePoint(ig));
157 }
158 764325 m_hasQuadPtCoor = true;
159 764325 }
160
161 /***********************************************************************************/
162 /***********************************************************************************/
163
164 19863341 double CurBaseFiniteElement::measure() const
165 {
166
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 19863341 times.
19863341 FEL_ASSERT(m_hasMeas);
167 19863341 double meas = 0.;
168 // The measure must not be computed from the boundary (keep m_numQuadraturePoint)!!
169
2/2
✓ Branch 1 taken 93492901 times.
✓ Branch 2 taken 19863341 times.
113356242 for ( int ig = 0; ig <m_numQuadraturePoint; ig++ ) meas +=weightMeas( ig );
170 19863341 return meas;
171 }
172
173 /***********************************************************************************/
174 /***********************************************************************************/
175
176 4319628 void CurBaseFiniteElement::coorMap(Point& pt, const Point& refpt) const
177 {
178 4319628 pt.clear();
179
2/2
✓ Branch 0 taken 12392808 times.
✓ Branch 1 taken 4319628 times.
16712436 for(int icoor=0; icoor<m_numCoor; icoor++) {
180
2/2
✓ Branch 0 taken 44145992 times.
✓ Branch 1 taken 12392808 times.
56538800 for (int i=0; i<m_numPoint; i++) {
181 44145992 pt[icoor] += m_point( i, icoor ) * m_geoEle.basisFunction().phi( i, refpt );
182 }
183 }
184 4319628 }
185
186 /***********************************************************************************/
187 /***********************************************************************************/
188
189 19588579 double CurBaseFiniteElement::diameter() const
190 {
191 19588579 double max = 0.;
192 19588579 double dist = 0.;
193
2/2
✓ Branch 0 taken 75149484 times.
✓ Branch 1 taken 19588579 times.
94738063 for (int i=0; i<m_numPoint; i++) {
194
2/2
✓ Branch 0 taken 107926053 times.
✓ Branch 1 taken 75149484 times.
183075537 for (int j=0; j<i; j++) {
195
2/2
✓ Branch 0 taken 321076770 times.
✓ Branch 1 taken 107926053 times.
429002823 for(int icoor=0; icoor<m_numCoor; icoor++) {
196 321076770 const double dp = m_point(i,icoor)-m_point(j,icoor);
197 321076770 dist += dp * dp;
198 }
199 107926053 dist = std::sqrt(dist);
200
2/2
✓ Branch 0 taken 42300925 times.
✓ Branch 1 taken 65625128 times.
107926053 if ( dist > max )
201 42300925 max = dist;
202 107926053 dist = 0.;
203 }
204 }
205 19588579 return max;
206 }
207
208 /***********************************************************************************/
209 /***********************************************************************************/
210
211 double CurBaseFiniteElement::measOfSegment(felInt id1, felInt id2) const
212 {
213 double dist=0;
214 for(int icoor=0; icoor<m_numCoor; icoor++) {
215 const double dp = m_point(id1,icoor)-m_point(id2,icoor);
216 dist += dp * dp;
217 }
218 dist = std::sqrt(dist);
219 return dist;
220 }
221
222 /***********************************************************************************/
223 /***********************************************************************************/
224
225 51628779 void CurBaseFiniteElement::updatePoint(const std::vector<Point*>& point)
226 {
227
2/2
✓ Branch 0 taken 196450304 times.
✓ Branch 1 taken 51628779 times.
248079083 for (int i=0; i<m_numPoint; i++) {
228
2/2
✓ Branch 0 taken 572260346 times.
✓ Branch 1 taken 196450304 times.
768710650 for(int icoor=0; icoor<m_numCoor; icoor++) {
229 572260346 m_point(i,icoor) = point[i]->coor(icoor);
230 }
231 }
232 51628779 m_hasPoint = true;
233 51628779 }
234
235 /***********************************************************************************/
236 /***********************************************************************************/
237
238 4080132 void CurBaseFiniteElement::updatePoint(const UBlasMatrix& point, const std::vector<int>& ipt)
239 {
240
2/2
✓ Branch 0 taken 8934012 times.
✓ Branch 1 taken 4080132 times.
13014144 for (int i=0; i<m_numPoint; i++) {
241
2/2
✓ Branch 0 taken 20189268 times.
✓ Branch 1 taken 8934012 times.
29123280 for(int icoor=0; icoor<m_numCoor; icoor++) {
242 20189268 m_point(i,icoor) = point(ipt[i],icoor);
243 }
244 }
245 4080132 m_hasPoint = true;
246 4080132 }
247
248 /***********************************************************************************/
249 /***********************************************************************************/
250
251 void CurBaseFiniteElement::update(const int id,const std::vector<Point*>& point)
252 {
253 updatePoint(point);
254 m_hasMeas = false;
255 m_hasQuadPtCoor = false;
256 m_currentId = id;
257 }
258
259 /***********************************************************************************/
260 /***********************************************************************************/
261
262 32888 void CurBaseFiniteElement::updateBasisWithQuadPoint(const std::vector<Point>& point)
263 {
264
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 32888 times.
32888 FEL_CHECK(m_numQuadraturePoint >= (felInt)point.size(), "Error when updating the basis function with new quadrature point: not enough points given");
265
266
2/2
✓ Branch 1 taken 80000 times.
✓ Branch 2 taken 32888 times.
112888 for(std::size_t ig=0; ig<point.size(); ig++) {
267
2/2
✓ Branch 0 taken 235888 times.
✓ Branch 1 taken 80000 times.
315888 for(int i=0; i<m_numDof; i++) {
268 235888 phi[ig](i) = m_refEle.basisFunction().phi(i, point[ig]);
269
2/2
✓ Branch 0 taken 463552 times.
✓ Branch 1 taken 235888 times.
699440 for(int icoor=0; icoor<m_numRefCoor; icoor++) {
270 463552 m_dPhiRef[ig](icoor, i) = m_refEle.basisFunction().dPhi(i, icoor, point[ig]);
271 }
272 }
273
2/2
✓ Branch 0 taken 235888 times.
✓ Branch 1 taken 80000 times.
315888 for(int k=0; k<m_numPoint; k++) {
274 235888 m_phiGeo[ig](k) = m_geoEle.basisFunction().phi(k, point[ig]);
275
2/2
✓ Branch 0 taken 463552 times.
✓ Branch 1 taken 235888 times.
699440 for(int icoor=0; icoor<m_numRefCoor; icoor++) {
276 463552 m_dPhiGeo[ig](icoor, k) = m_geoEle.basisFunction().dPhi(k, icoor, point[ig]);
277 }
278 }
279 }
280 32888 }
281
282 /***********************************************************************************/
283 /***********************************************************************************/
284
285 void CurBaseFiniteElement::mapPointInReferenceElement(const int id, const std::vector<Point*>& ptElem, const Point& inPoint, Point& outPoint)
286 {
287 updatePoint(ptElem);
288 m_hasMeas = false;
289 m_hasQuadPtCoor = false;
290 m_currentId = id;
291
292
293 // find the coordinate of the input points in the reference element of the element
294 switch(m_geoEle.shape().typeShape()) {
295 case Quadrilateral: {
296 // square
297 outPoint.x() = 2.*(inPoint.x() - ptElem[0]->x())/(ptElem[2]->x() - ptElem[0]->x()) - 1.;
298 outPoint.y() = 2.*(inPoint.y() - ptElem[0]->y())/(ptElem[2]->y() - ptElem[0]->y()) - 1.;
299 break;
300 }
301 case Triangle: {
302 // triangle
303 // this value can not be equal to zero unless the triangle is flat.
304 const double InvDet = 1./( (ptElem[2]->y() - ptElem[0]->y())*(ptElem[1]->x() - ptElem[0]->x())
305 - (ptElem[2]->x() - ptElem[0]->x())*(ptElem[1]->y() - ptElem[0]->y()) );
306
307 outPoint.x() = ( (ptElem[2]->y() - ptElem[0]->y()) * (inPoint.x() - ptElem[0]->x())
308 + (ptElem[0]->x() - ptElem[2]->x()) * (inPoint.y() - ptElem[0]->y()) ) * InvDet;
309
310 outPoint.y() = ( (ptElem[0]->y() - ptElem[1]->y()) * (inPoint.x() - ptElem[0]->x())
311 + (ptElem[1]->x() - ptElem[0]->x()) * (inPoint.y() - ptElem[0]->y()) ) * InvDet;
312 break;
313 }
314 case Tetrahedron: {
315 UBlasMatrix FLU(m_numRefCoor,m_numRefCoor);
316 UBlasMatrix inv_FLU(m_numRefCoor,m_numRefCoor);
317 UBlasVector refPoint(m_numRefCoor);
318 double det;
319
320 // Filling the matrix
321 for(int indx = 0; indx < m_numRefCoor ; ++indx){
322 FLU(0,indx) = ptElem[indx+1]->x() - ptElem[0]->x();
323 FLU(1,indx) = ptElem[indx+1]->y() - ptElem[0]->y();
324 FLU(2,indx) = ptElem[indx+1]->z() - ptElem[0]->z();
325 }
326
327 MathUtilities::InvertMatrix(FLU, inv_FLU, det);
328
329 // Filling the rhs
330 refPoint[0] = inPoint.x() - ptElem[0]->x();
331 refPoint[1] = inPoint.y() - ptElem[0]->y();
332 refPoint[2] = inPoint.z() - ptElem[0]->z();
333
334 refPoint = prod(inv_FLU,refPoint);
335
336 outPoint.x() = refPoint[0];
337 outPoint.y() = refPoint[1];
338 outPoint.z() = refPoint[2];
339
340 break;
341 }
342 case NullShape:
343 case Node:
344 case Segment:
345 case Hexahedron:
346 case Prism:
347 case Pyramid: {
348 FEL_ERROR("The transformation is not implemented for this type of mesh cells");
349 break;
350 }
351 }
352 }
353
354 /***********************************************************************************/
355 /***********************************************************************************/
356
357 void CurBaseFiniteElement::m_printBase(int verbose,std::ostream& c) const
358 {
359 if(verbose) {
360 c << " RefElement: " << m_refEle.name() << "\n";
361 c << " GeoElement: " << m_geoEle.name() << "\n";
362 c << " QuadratureRule: " << m_quadratureRule->name() << "\n";
363 c << " numDof: " << m_numDof << "\n";
364 c << " numQuadPoint: " << m_numQuadraturePoint << "\n";
365 c << " numRefCoor: " << m_numRefCoor << "\n";
366 c << " numCoor: " << m_numCoor << "\n";
367 c << " numPoint: " << m_numPoint << "\n";
368 c << " currentID: " << m_currentId << "\n";
369 /*
370 if( m_hasPoint && verbose > 1 ) c << " point = " << m_point << "\n";
371 if( m_hasMeas && verbose > 1 ) c << " meas = " << measure() << "\n";
372 if( verbose > 2 ) {
373 if ( theNumQuadPointTotal != m_numQuadraturePoint ) {
374 c << " The numQuadPoint=" << m_numQuadraturePoint << " first are internal. The others live on the boundary elements of the element.\n";}
375 for (int ig = 0;ig <theNumQuadPointTotal;ig++ ){
376 c << " phi[ig=" <<ig<< "] = " << phi[ig] << "\n";
377 }
378 }
379 if( verbose > 2 ) {
380 for (int ig = 0;ig <theNumQuadPointTotal;ig++ ){
381 c << " phiGeo[ig=" <<ig<< "] = " << m_phiGeo[ig] << "\n";
382 }
383 }
384 if( m_hasFirstDeriv && verbose > 2 ) {
385 for (int ig = 0;ig <theNumQuadPointTotal;ig++ ){
386 c << " dphi[ig=" <<ig<< "] = " << dPhi[ig] << "\n";
387 }
388 }
389 if( verbose > 2 ) {
390 for (int ig = 0;ig <theNumQuadPointTotal;ig++ ){
391 c << " dphiGeo[ig=" <<ig<< "] = " << m_dPhiGeo[ig] << "\n";
392 }
393 }
394 if( verbose > 2 ) {
395 for (int ig = 0;ig <theNumQuadPointTotal;ig++ ){
396 c << " dphiRef[ig=" <<ig<< "] = " << m_dPhiRef[ig] << "\n";
397 }
398 }
399 if( m_hasFirstDeriv && verbose > 2 ) {
400 for (int ig = 0;ig <theNumQuadPointTotal;ig++ ){
401 c << " jacobian[ig=" <<ig<< "] = " << m_jacobian[ig] << "\n";
402 }
403 }
404 if( m_hasQuadPtCoor && verbose > 2 ) {
405 for (int ig = 0;ig <theNumQuadPointTotal;ig++ ){
406 c << " currentQuadPoint[ig=" <<ig<< "] : "; currentQuadPoint[ig].print(verbose,c);
407 }
408 }
409 */
410 }
411 }
412
413 /***********************************************************************************/
414 /***********************************************************************************/
415
416 void CurBaseFiniteElement::print(int verbose,std::ostream& c) const
417 {
418 if(verbose) {
419 c << "CurBaseFiniteElement: " << "\n";
420 m_printBase(verbose, c);
421 std::cout << std::endl;
422 }
423 }
424 }
425