Dune Core Modules (2.5.0)

affinegeometry.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_GEOMETRY_AFFINEGEOMETRY_HH
4#define DUNE_GEOMETRY_AFFINEGEOMETRY_HH
5
11#include <cmath>
12
15
16#include <dune/geometry/type.hh>
17
18namespace Dune
19{
20
21 // External Forward Declarations
22 // -----------------------------
23
24 template< class ctype, int dim >
25 class ReferenceElement;
26
27 template< class ctype, int dim >
28 struct ReferenceElements;
29
30
31
32 namespace Impl
33 {
34
35 // FieldMatrixHelper
36 // -----------------
37
38 template< class ct >
39 struct FieldMatrixHelper
40 {
41 typedef ct ctype;
42
43 template< int m, int n >
44 static void Ax ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, n > &x, FieldVector< ctype, m > &ret )
45 {
46 for( int i = 0; i < m; ++i )
47 {
48 ret[ i ] = ctype( 0 );
49 for( int j = 0; j < n; ++j )
50 ret[ i ] += A[ i ][ j ] * x[ j ];
51 }
52 }
53
54 template< int m, int n >
55 static void ATx ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, m > &x, FieldVector< ctype, n > &ret )
56 {
57 for( int i = 0; i < n; ++i )
58 {
59 ret[ i ] = ctype( 0 );
60 for( int j = 0; j < m; ++j )
61 ret[ i ] += A[ j ][ i ] * x[ j ];
62 }
63 }
64
65 template< int m, int n, int p >
66 static void AB ( const FieldMatrix< ctype, m, n > &A, const FieldMatrix< ctype, n, p > &B, FieldMatrix< ctype, m, p > &ret )
67 {
68 for( int i = 0; i < m; ++i )
69 {
70 for( int j = 0; j < p; ++j )
71 {
72 ret[ i ][ j ] = ctype( 0 );
73 for( int k = 0; k < n; ++k )
74 ret[ i ][ j ] += A[ i ][ k ] * B[ k ][ j ];
75 }
76 }
77 }
78
79 template< int m, int n, int p >
80 static void ATBT ( const FieldMatrix< ctype, m, n > &A, const FieldMatrix< ctype, p, m > &B, FieldMatrix< ctype, n, p > &ret )
81 {
82 for( int i = 0; i < n; ++i )
83 {
84 for( int j = 0; j < p; ++j )
85 {
86 ret[ i ][ j ] = ctype( 0 );
87 for( int k = 0; k < m; ++k )
88 ret[ i ][ j ] += A[ k ][ i ] * B[ j ][ k ];
89 }
90 }
91 }
92
93 template< int m, int n >
94 static void ATA_L ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, n, n > &ret )
95 {
96 for( int i = 0; i < n; ++i )
97 {
98 for( int j = 0; j <= i; ++j )
99 {
100 ret[ i ][ j ] = ctype( 0 );
101 for( int k = 0; k < m; ++k )
102 ret[ i ][ j ] += A[ k ][ i ] * A[ k ][ j ];
103 }
104 }
105 }
106
107 template< int m, int n >
108 static void ATA ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, n, n > &ret )
109 {
110 for( int i = 0; i < n; ++i )
111 {
112 for( int j = 0; j <= i; ++j )
113 {
114 ret[ i ][ j ] = ctype( 0 );
115 for( int k = 0; k < m; ++k )
116 ret[ i ][ j ] += A[ k ][ i ] * A[ k ][ j ];
117 ret[ j ][ i ] = ret[ i ][ j ];
118 }
119
120 ret[ i ][ i ] = ctype( 0 );
121 for( int k = 0; k < m; ++k )
122 ret[ i ][ i ] += A[ k ][ i ] * A[ k ][ i ];
123 }
124 }
125
126 template< int m, int n >
127 static void AAT_L ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, m, m > &ret )
128 {
129 /*
130 if (m==2) {
131 ret[0][0] = A[0]*A[0];
132 ret[1][1] = A[1]*A[1];
133 ret[1][0] = A[0]*A[1];
134 }
135 else
136 */
137 for( int i = 0; i < m; ++i )
138 {
139 for( int j = 0; j <= i; ++j )
140 {
141 ctype &retij = ret[ i ][ j ];
142 retij = A[ i ][ 0 ] * A[ j ][ 0 ];
143 for( int k = 1; k < n; ++k )
144 retij += A[ i ][ k ] * A[ j ][ k ];
145 }
146 }
147 }
148
149 template< int m, int n >
150 static void AAT ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, m, m > &ret )
151 {
152 for( int i = 0; i < m; ++i )
153 {
154 for( int j = 0; j < i; ++j )
155 {
156 ret[ i ][ j ] = ctype( 0 );
157 for( int k = 0; k < n; ++k )
158 ret[ i ][ j ] += A[ i ][ k ] * A[ j ][ k ];
159 ret[ j ][ i ] = ret[ i ][ j ];
160 }
161 ret[ i ][ i ] = ctype( 0 );
162 for( int k = 0; k < n; ++k )
163 ret[ i ][ i ] += A[ i ][ k ] * A[ i ][ k ];
164 }
165 }
166
167 template< int n >
168 static void Lx ( const FieldMatrix< ctype, n, n > &L, const FieldVector< ctype, n > &x, FieldVector< ctype, n > &ret )
169 {
170 for( int i = 0; i < n; ++i )
171 {
172 ret[ i ] = ctype( 0 );
173 for( int j = 0; j <= i; ++j )
174 ret[ i ] += L[ i ][ j ] * x[ j ];
175 }
176 }
177
178 template< int n >
179 static void LTx ( const FieldMatrix< ctype, n, n > &L, const FieldVector< ctype, n > &x, FieldVector< ctype, n > &ret )
180 {
181 for( int i = 0; i < n; ++i )
182 {
183 ret[ i ] = ctype( 0 );
184 for( int j = i; j < n; ++j )
185 ret[ i ] += L[ j ][ i ] * x[ j ];
186 }
187 }
188
189 template< int n >
190 static void LTL ( const FieldMatrix< ctype, n, n > &L, FieldMatrix< ctype, n, n > &ret )
191 {
192 for( int i = 0; i < n; ++i )
193 {
194 for( int j = 0; j < i; ++j )
195 {
196 ret[ i ][ j ] = ctype( 0 );
197 for( int k = i; k < n; ++k )
198 ret[ i ][ j ] += L[ k ][ i ] * L[ k ][ j ];
199 ret[ j ][ i ] = ret[ i ][ j ];
200 }
201 ret[ i ][ i ] = ctype( 0 );
202 for( int k = i; k < n; ++k )
203 ret[ i ][ i ] += L[ k ][ i ] * L[ k ][ i ];
204 }
205 }
206
207 template< int n >
208 static void LLT ( const FieldMatrix< ctype, n, n > &L, FieldMatrix< ctype, n, n > &ret )
209 {
210 for( int i = 0; i < n; ++i )
211 {
212 for( int j = 0; j < i; ++j )
213 {
214 ret[ i ][ j ] = ctype( 0 );
215 for( int k = 0; k <= j; ++k )
216 ret[ i ][ j ] += L[ i ][ k ] * L[ j ][ k ];
217 ret[ j ][ i ] = ret[ i ][ j ];
218 }
219 ret[ i ][ i ] = ctype( 0 );
220 for( int k = 0; k <= i; ++k )
221 ret[ i ][ i ] += L[ i ][ k ] * L[ i ][ k ];
222 }
223 }
224
225 template< int n >
226 static void cholesky_L ( const FieldMatrix< ctype, n, n > &A, FieldMatrix< ctype, n, n > &ret )
227 {
228 for( int i = 0; i < n; ++i )
229 {
230 ctype &rii = ret[ i ][ i ];
231
232 ctype xDiag = A[ i ][ i ];
233 for( int j = 0; j < i; ++j )
234 xDiag -= ret[ i ][ j ] * ret[ i ][ j ];
235 assert( xDiag > ctype( 0 ) );
236 rii = sqrt( xDiag );
237
238 ctype invrii = ctype( 1 ) / rii;
239 for( int k = i+1; k < n; ++k )
240 {
241 ctype x = A[ k ][ i ];
242 for( int j = 0; j < i; ++j )
243 x -= ret[ i ][ j ] * ret[ k ][ j ];
244 ret[ k ][ i ] = invrii * x;
245 }
246 }
247 }
248
249 template< int n >
250 static ctype detL ( const FieldMatrix< ctype, n, n > &L )
251 {
252 ctype det( 1 );
253 for( int i = 0; i < n; ++i )
254 det *= L[ i ][ i ];
255 return det;
256 }
257
258 template< int n >
259 static ctype invL ( FieldMatrix< ctype, n, n > &L )
260 {
261 ctype det( 1 );
262 for( int i = 0; i < n; ++i )
263 {
264 ctype &lii = L[ i ][ i ];
265 det *= lii;
266 lii = ctype( 1 ) / lii;
267 for( int j = 0; j < i; ++j )
268 {
269 ctype &lij = L[ i ][ j ];
270 ctype x = lij * L[ j ][ j ];
271 for( int k = j+1; k < i; ++k )
272 x += L[ i ][ k ] * L[ k ][ j ];
273 lij = (-lii) * x;
274 }
275 }
276 return det;
277 }
278
279 // calculates x := L^{-1} x
280 template< int n >
281 static void invLx ( FieldMatrix< ctype, n, n > &L, FieldVector< ctype, n > &x )
282 {
283 for( int i = 0; i < n; ++i )
284 {
285 for( int j = 0; j < i; ++j )
286 x[ i ] -= L[ i ][ j ] * x[ j ];
287 x[ i ] /= L[ i ][ i ];
288 }
289 }
290
291 // calculates x := L^{-T} x
292 template< int n >
293 static void invLTx ( FieldMatrix< ctype, n, n > &L, FieldVector< ctype, n > &x )
294 {
295 for( int i = n; i > 0; --i )
296 {
297 for( int j = i; j < n; ++j )
298 x[ i-1 ] -= L[ j ][ i-1 ] * x[ j ];
299 x[ i-1 ] /= L[ i-1 ][ i-1 ];
300 }
301 }
302
303 template< int n >
304 static ctype spdDetA ( const FieldMatrix< ctype, n, n > &A )
305 {
306 // return A[0][0]*A[1][1]-A[1][0]*A[1][0];
307 FieldMatrix< ctype, n, n > L;
308 cholesky_L( A, L );
309 return detL( L );
310 }
311
312 template< int n >
313 static ctype spdInvA ( FieldMatrix< ctype, n, n > &A )
314 {
315 FieldMatrix< ctype, n, n > L;
316 cholesky_L( A, L );
317 const ctype det = invL( L );
318 LTL( L, A );
319 return det;
320 }
321
322 // calculate x := A^{-1} x
323 template< int n >
324 static void spdInvAx ( FieldMatrix< ctype, n, n > &A, FieldVector< ctype, n > &x )
325 {
326 FieldMatrix< ctype, n, n > L;
327 cholesky_L( A, L );
328 invLx( L, x );
329 invLTx( L, x );
330 }
331
332 template< int m, int n >
333 static ctype detATA ( const FieldMatrix< ctype, m, n > &A )
334 {
335 if( m >= n )
336 {
337 FieldMatrix< ctype, n, n > ata;
338 ATA_L( A, ata );
339 return spdDetA( ata );
340 }
341 else
342 return ctype( 0 );
343 }
344
350 template< int m, int n >
351 static ctype sqrtDetAAT ( const FieldMatrix< ctype, m, n > &A )
352 {
353 using std::abs;
354 using std::sqrt;
355 // These special cases are here not only for speed reasons:
356 // The general implementation aborts if the matrix is almost singular,
357 // and the special implementation provide a stable way to handle that case.
358 if( (n == 2) && (m == 2) )
359 {
360 // Special implementation for 2x2 matrices: faster and more stable
361 return abs( A[ 0 ][ 0 ]*A[ 1 ][ 1 ] - A[ 1 ][ 0 ]*A[ 0 ][ 1 ] );
362 }
363 else if( (n == 3) && (m == 3) )
364 {
365 // Special implementation for 3x3 matrices
366 const ctype v0 = A[ 0 ][ 1 ] * A[ 1 ][ 2 ] - A[ 1 ][ 1 ] * A[ 0 ][ 2 ];
367 const ctype v1 = A[ 0 ][ 2 ] * A[ 1 ][ 0 ] - A[ 1 ][ 2 ] * A[ 0 ][ 0 ];
368 const ctype v2 = A[ 0 ][ 0 ] * A[ 1 ][ 1 ] - A[ 1 ][ 0 ] * A[ 0 ][ 1 ];
369 return abs( v0 * A[ 2 ][ 0 ] + v1 * A[ 2 ][ 1 ] + v2 * A[ 2 ][ 2 ] );
370 }
371 else if ( (n == 3) && (m == 2) )
372 {
373 // Special implementation for 2x3 matrices
374 const ctype v0 = A[ 0 ][ 0 ] * A[ 1 ][ 1 ] - A[ 0 ][ 1 ] * A[ 1 ][ 0 ];
375 const ctype v1 = A[ 0 ][ 0 ] * A[ 1 ][ 2 ] - A[ 1 ][ 0 ] * A[ 0 ][ 2 ];
376 const ctype v2 = A[ 0 ][ 1 ] * A[ 1 ][ 2 ] - A[ 0 ][ 2 ] * A[ 1 ][ 1 ];
377 return sqrt( v0*v0 + v1*v1 + v2*v2);
378 }
379 else if( n >= m )
380 {
381 // General case
382 FieldMatrix< ctype, m, m > aat;
383 AAT_L( A, aat );
384 return spdDetA( aat );
385 }
386 else
387 return ctype( 0 );
388 }
389
390 // A^{-1}_L = (A^T A)^{-1} A^T
391 // => A^{-1}_L A = I
392 template< int m, int n >
393 static ctype leftInvA ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, n, m > &ret )
394 {
395 static_assert((m >= n), "Matrix has no left inverse.");
396 FieldMatrix< ctype, n, n > ata;
397 ATA_L( A, ata );
398 const ctype det = spdInvA( ata );
399 ATBT( ata, A, ret );
400 return det;
401 }
402
403 template< int m, int n >
404 static void leftInvAx ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, m > &x, FieldVector< ctype, n > &y )
405 {
406 static_assert((m >= n), "Matrix has no left inverse.");
407 FieldMatrix< ctype, n, n > ata;
408 ATx( A, x, y );
409 ATA_L( A, ata );
410 spdInvAx( ata, y );
411 }
412
414 template< int m, int n >
415 static ctype rightInvA ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, n, m > &ret )
416 {
417 static_assert((n >= m), "Matrix has no right inverse.");
418 using std::abs;
419 if( (n == 2) && (m == 2) )
420 {
421 const ctype det = (A[ 0 ][ 0 ]*A[ 1 ][ 1 ] - A[ 1 ][ 0 ]*A[ 0 ][ 1 ]);
422 const ctype detInv = ctype( 1 ) / det;
423 ret[ 0 ][ 0 ] = A[ 1 ][ 1 ] * detInv;
424 ret[ 1 ][ 1 ] = A[ 0 ][ 0 ] * detInv;
425 ret[ 1 ][ 0 ] = -A[ 1 ][ 0 ] * detInv;
426 ret[ 0 ][ 1 ] = -A[ 0 ][ 1 ] * detInv;
427 return abs( det );
428 }
429 else
430 {
431 FieldMatrix< ctype, m , m > aat;
432 AAT_L( A, aat );
433 const ctype det = spdInvA( aat );
434 ATBT( A , aat , ret );
435 return det;
436 }
437 }
438
439 template< int m, int n >
440 static void xTRightInvA ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, n > &x, FieldVector< ctype, m > &y )
441 {
442 static_assert((n >= m), "Matrix has no right inverse.");
443 FieldMatrix< ctype, m, m > aat;
444 Ax( A, x, y );
445 AAT_L( A, aat );
446 spdInvAx( aat, y );
447 }
448 };
449
450 } // namespace Impl
451
452
453
459 template< class ct, int mydim, int cdim>
461 {
462 public:
463
465 typedef ct ctype;
466
468 static const int mydimension= mydim;
469
471 static const int coorddimension = cdim;
472
475
478
481
484
485 private:
488
490
491 // Helper class to compute a matrix pseudo inverse
492 typedef Impl::FieldMatrixHelper< ct > MatrixHelper;
493
494 public:
496 AffineGeometry ( const ReferenceElement &refElement, const GlobalCoordinate &origin,
497 const JacobianTransposed &jt )
498 : refElement_(&refElement), origin_(origin), jacobianTransposed_(jt)
499 {
500 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
501 }
502
505 const JacobianTransposed &jt )
506 : refElement_( &ReferenceElements::general( gt ) ), origin_(origin), jacobianTransposed_( jt )
507 {
508 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
509 }
510
512 template< class CoordVector >
513 AffineGeometry ( const ReferenceElement &refElement, const CoordVector &coordVector )
514 : refElement_(&refElement), origin_(coordVector[0])
515 {
516 for( int i = 0; i < mydimension; ++i )
517 jacobianTransposed_[ i ] = coordVector[ i+1 ] - origin_;
518 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
519 }
520
522 template< class CoordVector >
523 AffineGeometry ( Dune::GeometryType gt, const CoordVector &coordVector )
524 : refElement_(&ReferenceElements::general( gt )), origin_(coordVector[0] )
525 {
526 for( int i = 0; i < mydimension; ++i )
527 jacobianTransposed_[ i ] = coordVector[ i+1 ] - origin_;
528 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
529 }
530
532 bool affine () const { return true; }
533
535 Dune::GeometryType type () const { return refElement_->type(); }
536
538 int corners () const { return refElement_->size( mydimension ); }
539
541 GlobalCoordinate corner ( int i ) const
542 {
543 return global( refElement_->position( i, mydimension ) );
544 }
545
547 GlobalCoordinate center () const { return global( refElement_->position( 0, 0 ) ); }
548
556 {
557 GlobalCoordinate global( origin_ );
558 jacobianTransposed_.umtv( local, global );
559 return global;
560 }
561
576 {
577 LocalCoordinate local;
578 jacobianInverseTransposed_.mtv( global - origin_, local );
579 return local;
580 }
581
593 {
595 return integrationElement_;
596 }
597
599 ctype volume () const
600 {
601 return integrationElement_ * refElement_->volume();
602 }
603
611 {
613 return jacobianTransposed_;
614 }
615
623 {
625 return jacobianInverseTransposed_;
626 }
627
628 friend const ReferenceElement &referenceElement ( const AffineGeometry &geometry ) { return *geometry.refElement_; }
629
630 private:
631 const ReferenceElement* refElement_;
632 GlobalCoordinate origin_;
633 JacobianTransposed jacobianTransposed_;
634 JacobianInverseTransposed jacobianInverseTransposed_;
635 ctype integrationElement_;
636 };
637
638} // namespace Dune
639
640#endif // #ifndef DUNE_GEOMETRY_AFFINEGEOMETRY_HH
Implementation of the Geometry interface for affine geometries.
Definition: affinegeometry.hh:461
AffineGeometry(const ReferenceElement &refElement, const CoordVector &coordVector)
Create affine geometry from reference element and a vector of vertex coordinates.
Definition: affinegeometry.hh:513
AffineGeometry(Dune::GeometryType gt, const GlobalCoordinate &origin, const JacobianTransposed &jt)
Create affine geometry from GeometryType, one vertex, and the Jacobian matrix.
Definition: affinegeometry.hh:504
FieldVector< ctype, mydimension > LocalCoordinate
Type for local coordinate vector.
Definition: affinegeometry.hh:474
Dune::GeometryType type() const
Obtain the type of the reference element.
Definition: affinegeometry.hh:535
static const int mydimension
Dimension of the geometry.
Definition: affinegeometry.hh:468
AffineGeometry(const ReferenceElement &refElement, const GlobalCoordinate &origin, const JacobianTransposed &jt)
Create affine geometry from reference element, one vertex, and the Jacobian matrix.
Definition: affinegeometry.hh:496
ctype volume() const
Obtain the volume of the element.
Definition: affinegeometry.hh:599
AffineGeometry(Dune::GeometryType gt, const CoordVector &coordVector)
Create affine geometry from GeometryType and a vector of vertex coordinates.
Definition: affinegeometry.hh:523
ctype integrationElement(const LocalCoordinate &local) const
Obtain the integration element.
Definition: affinegeometry.hh:592
const JacobianInverseTransposed & jacobianInverseTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian's inverse.
Definition: affinegeometry.hh:622
FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
Type for the transposed Jacobian matrix.
Definition: affinegeometry.hh:480
GlobalCoordinate corner(int i) const
Obtain coordinates of the i-th corner.
Definition: affinegeometry.hh:541
int corners() const
Obtain number of corners of the corresponding reference element.
Definition: affinegeometry.hh:538
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed
Type for the transposed inverse Jacobian matrix.
Definition: affinegeometry.hh:483
static const int coorddimension
Dimension of the world space.
Definition: affinegeometry.hh:471
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the mapping.
Definition: affinegeometry.hh:555
GlobalCoordinate center() const
Obtain the centroid of the mapping's image.
Definition: affinegeometry.hh:547
ct ctype
Type used for coordinates.
Definition: affinegeometry.hh:465
FieldVector< ctype, coorddimension > GlobalCoordinate
Type for coordinate vector in world space.
Definition: affinegeometry.hh:477
bool affine() const
Always true: this is an affine geometry.
Definition: affinegeometry.hh:532
const JacobianTransposed & jacobianTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian.
Definition: affinegeometry.hh:610
void mtv(const X &x, Y &y) const
y = A^T x
Definition: densematrix.hh:380
void umtv(const X &x, Y &y) const
y += A^T x
Definition: densematrix.hh:408
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:268
This class provides access to geometric and topological properties of a reference element.
Definition: referenceelements.hh:354
int size(int c) const
number of subentities of codimension c
Definition: referenceelements.hh:381
const FieldVector< ctype, dim > & position(int i, int c) const
position of the barycenter of entity (i,c)
Definition: referenceelements.hh:449
ctype volume() const
obtain the volume of the reference element
Definition: referenceelements.hh:486
const GeometryType & type(int i, int c) const
obtain the type of subentity (i,c)
Definition: referenceelements.hh:431
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:132
Dune namespace.
Definition: alignment.hh:11
Class providing access to the singletons of the reference elements.
Definition: referenceelements.hh:752
A unique label for each type of element that can occur in a grid.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)