OpenWalnut  1.5.0dev
WDataSetSphericalHarmonics.cpp
1 //---------------------------------------------------------------------------
2 //
3 // Project: OpenWalnut ( http://www.openwalnut.org )
4 //
5 // Copyright 2009 OpenWalnut Community, BSV@Uni-Leipzig and CNCF@MPI-CBS
6 // For more information see http://www.openwalnut.org/copying
7 //
8 // This file is part of OpenWalnut.
9 //
10 // OpenWalnut is free software: you can redistribute it and/or modify
11 // it under the terms of the GNU Lesser General Public License as published by
12 // the Free Software Foundation, either version 3 of the License, or
13 // (at your option) any later version.
14 //
15 // OpenWalnut is distributed in the hope that it will be useful,
16 // but WITHOUT ANY WARRANTY; without even the implied warranty of
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 // GNU Lesser General Public License for more details.
19 //
20 // You should have received a copy of the GNU Lesser General Public License
21 // along with OpenWalnut. If not, see <http://www.gnu.org/licenses/>.
22 //
23 //---------------------------------------------------------------------------
24 
25 #include <cmath>
26 #include <memory>
27 #include <stdint.h>
28 #include <string>
29 #include <vector>
30 
31 #include "../common/WAssert.h"
32 #include "../common/math/WSymmetricSphericalHarmonic.h"
33 #include "../common/math/linearAlgebra/WPosition.h"
34 #include "WDataSetSingle.h"
35 #include "WDataSetSphericalHarmonics.h"
36 
37 // prototype instance as singleton
38 std::shared_ptr< WPrototyped > WDataSetSphericalHarmonics::m_prototype = std::shared_ptr< WPrototyped >();
39 
40 WDataSetSphericalHarmonics::WDataSetSphericalHarmonics( std::shared_ptr< WValueSetBase > newValueSet,
41  std::shared_ptr< WGrid > newGrid ) :
42  WDataSetSingle( newValueSet, newGrid ), m_valueSet( newValueSet )
43 {
44  m_gridRegular3D = std::dynamic_pointer_cast< WGridRegular3D >( newGrid );
45  WAssert( newValueSet, "No value set given." );
46  WAssert( newGrid, "No grid given." );
47 }
48 
50  : WDataSetSingle()
51 {
52 }
53 
55 {
56 }
57 
58 WDataSetSingle::SPtr WDataSetSphericalHarmonics::clone( std::shared_ptr< WValueSetBase > newValueSet, std::shared_ptr< WGrid > newGrid ) const
59 {
60  return WDataSetSingle::SPtr( new WDataSetSphericalHarmonics( newValueSet, newGrid ) );
61 }
62 
63 WDataSetSingle::SPtr WDataSetSphericalHarmonics::clone( std::shared_ptr< WValueSetBase > newValueSet ) const
64 {
65  return WDataSetSingle::SPtr( new WDataSetSphericalHarmonics( newValueSet, getGrid() ) );
66 }
67 
68 WDataSetSingle::SPtr WDataSetSphericalHarmonics::clone( std::shared_ptr< WGrid > newGrid ) const
69 {
71 }
72 
74 {
76 }
77 
78 std::shared_ptr< WPrototyped > WDataSetSphericalHarmonics::getPrototype()
79 {
80  if( !m_prototype )
81  {
82  m_prototype = std::shared_ptr< WPrototyped >( new WDataSetSphericalHarmonics() );
83  }
84 
85  return m_prototype;
86 }
87 
89 {
90  *success = m_gridRegular3D->encloses( pos );
91 
92  bool isInside = true;
93  size_t cellId = m_gridRegular3D->getCellId( pos, &isInside );
94 
95  if( !isInside )
96  {
97  *success = false;
99  }
100 
101  // ids of vertices for interpolation
102  WGridRegular3D::CellVertexArray vertexIds = m_gridRegular3D->getCellVertexIds( cellId );
103 
104  WPosition localPos = m_gridRegular3D->getTransform().directionToGridSpace( pos - m_gridRegular3D->getPosition( vertexIds[0] ) );
105 
106  double lambdaX = localPos[0];
107  double lambdaY = localPos[1];
108  double lambdaZ = localPos[2];
109 
110  WValue< double > h( 8 );
111 // lZ lY
112 // | /
113 // | 6___/_7
114 // |/: /|
115 // 4_:___5 |
116 // | :...|.|
117 // |.2 | 3
118 // |_____|/ ____lX
119 // 0 1
120  h[0] = ( 1 - lambdaX ) * ( 1 - lambdaY ) * ( 1 - lambdaZ );
121  h[1] = ( lambdaX ) * ( 1 - lambdaY ) * ( 1 - lambdaZ );
122  h[2] = ( 1 - lambdaX ) * ( lambdaY ) * ( 1 - lambdaZ );
123  h[3] = ( lambdaX ) * ( lambdaY ) * ( 1 - lambdaZ );
124  h[4] = ( 1 - lambdaX ) * ( 1 - lambdaY ) * ( lambdaZ );
125  h[5] = ( lambdaX ) * ( 1 - lambdaY ) * ( lambdaZ );
126  h[6] = ( 1 - lambdaX ) * ( lambdaY ) * ( lambdaZ );
127  h[7] = ( lambdaX ) * ( lambdaY ) * ( lambdaZ );
128 
129  // take
130  WValue<double> interpolatedCoefficients( m_valueSet->dimension() );
131  for( size_t i = 0; i < 8; ++i )
132  {
133  interpolatedCoefficients += h[i] * m_valueSet->getWValueDouble( vertexIds[i] );
134  }
135 
136  *success = true;
137 
138  return WSymmetricSphericalHarmonic< double >( interpolatedCoefficients );
139 }
140 
142 {
143  if( index < m_valueSet->size() ) return WSymmetricSphericalHarmonic< double >( m_valueSet->getWValueDouble( index ) );
144  wlog::warn( "WDataSetSphericalHarmonics" ) << "Invalid index. Must be smaller than " << m_valueSet->size();
146 }
147 
148 const std::string WDataSetSphericalHarmonics::getName() const
149 {
150  return "WDataSetSphericalHarmonics";
151 }
152 
154 {
155  return "Contains factors for spherical harmonics.";
156 }
157 
159 {
160  return false;
161 }
162 
A data set consisting of a set of values based on a grid.
std::shared_ptr< WValueSetBase > getValueSet() const
std::shared_ptr< WGrid > getGrid() const
std::shared_ptr< WDataSetSingle > SPtr
Convenience typedef for a std::shared_ptr.
std::shared_ptr< WValueSetBase > m_valueSet
The valueset of the data set.
virtual const std::string getDescription() const
Gets the description for this prototype.
WSymmetricSphericalHarmonic< double > interpolate(const WPosition &pos, bool *success) const
Interpolates the field of spherical harmonics at the given position.
virtual ~WDataSetSphericalHarmonics()
Destroys this DataSet instance.
static std::shared_ptr< WPrototyped > m_prototype
The prototype as singleton.
virtual WDataSetSingle::SPtr clone() const
Creates a copy (clone) of this instance.
virtual const std::string getName() const
Gets the name of this prototype.
std::shared_ptr< WGridRegular3D > m_gridRegular3D
The regular 3d grid of the data set.
static std::shared_ptr< WPrototyped > getPrototype()
Returns a prototype instantiated with the true type of the deriving class.
virtual bool isTexture() const
Determines whether this dataset can be used as a texture.
WDataSetSphericalHarmonics()
Construct an empty and unusable instance.
WSymmetricSphericalHarmonic< double > getSphericalHarmonicAt(size_t index) const
Get the spherical harmonic on the given position in value set.
boost::array< size_t, 8 > CellVertexArray
Convenience typedef for a boost::array< size_t, 8 >.
This only is a 3d double vector.
Class for symmetric spherical harmonics The index scheme of the coefficients/basis values is like in ...
Base class for all higher level values like tensors, vectors, matrices and so on.
Definition: WValue.h:41
WStreamedLogger warn(const std::string &source)
Logging a warning message.
Definition: WLogger.h:309