InitContactsForHCSITS.h 7.1 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
//======================================================================================================================
//
//  This file is part of waLBerla. waLBerla is free software: you can
//  redistribute it and/or modify it under the terms of the GNU General Public
//  License as published by the Free Software Foundation, either version 3 of
//  the License, or (at your option) any later version.
//
//  waLBerla is distributed in the hope that it will be useful, but WITHOUT
//  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
//  FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
//  for more details.
//
//  You should have received a copy of the GNU General Public License along
//  with waLBerla (see COPYING.txt). If not, see <http://www.gnu.org/licenses/>.
//
16
//! \file
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
//! \author Tobias Leemann <tobias.leemann@fau.de>
//
//======================================================================================================================

//======================================================================================================================
//
//  THIS FILE IS GENERATED - PLEASE CHANGE THE TEMPLATE !!!
//
//======================================================================================================================

#pragma once

#include <mesa_pd/common/ParticleFunctions.h>
#include <mesa_pd/data/DataTypes.h>
#include <mesa_pd/data/IAccessor.h>
#include <mesa_pd/data/IContactAccessor.h>
#include <mesa_pd/data/Flags.h>
#include <core/logging/Logging.h>

#include <vector>

namespace walberla {
namespace mesa_pd {
namespace kernel {

/**
* Init the datastructures for a contact for later use of the HCSITS-Solver.
* Call this kernel on all contacts that you want to treat with HCSITS before performing any relaxation timesteps.
* Use setErp() to set the error reduction parameter, which defines the share of the overlap that is resolved in one step
* (0 meaning after the relaxation the overlap is the same, 1 meaning the particles will have no overlap after relaxation).
* Use setFriction(a,b, cof) to define the coefficient of friction cof between materials a and b. It is assumed to be
 * symmetric w.r.t. the materials.
* \ingroup mesa_pd_kernel
*/
class InitContactsForHCSITS{
public:

   // Default constructor sets the default values
   explicit InitContactsForHCSITS(const uint_t numParticleTypes) :
   numParticleTypes_ (numParticleTypes),
   erp_( real_t(0.8) ),
   maximumPenetration_( 0 )
   {

      
      friction_.resize(numParticleTypes * numParticleTypes, real_t(0));
   }

   // Getter and Setter Functions
   inline const real_t& getErp() const {return erp_;}
   inline void setErp(real_t v) { erp_ = v;}
   
   inline const real_t& getMaximumPenetration() const {return maximumPenetration_;}
   inline void setMaximumPenetration(real_t v) { maximumPenetration_ = v;}
   

   // Getter and Setter Functions for material parameter combinations (they are assumed symmetric).
   
   inline void setFriction(const size_t type1, const size_t type2, const real_t& val)
   {
      WALBERLA_ASSERT_LESS( type1, numParticleTypes_ );
      WALBERLA_ASSERT_LESS( type2, numParticleTypes_ );
      friction_[numParticleTypes_*type1 + type2] = val;
      friction_[numParticleTypes_*type2 + type1] = val;
   }

   
   inline real_t getFriction(const size_t type1, const size_t type2) const
   {
      WALBERLA_ASSERT_LESS( type1, numParticleTypes_ );
      WALBERLA_ASSERT_LESS( type2, numParticleTypes_ );
      WALBERLA_ASSERT_FLOAT_EQUAL( friction_[numParticleTypes_*type1 + type2],
              friction_[numParticleTypes_*type2 + type1],
              "parameter matrix for friction not symmetric!");
      return friction_[numParticleTypes_*type1 + type2];
   }

   // List material parameters
   
   std::vector<real_t> friction_ {};

   template<typename CAccessor, typename PAccessor>
   void operator()(size_t c, CAccessor &ca, PAccessor &pa);

private:
   // List properties
   const uint_t numParticleTypes_;

   
   real_t erp_;
   real_t maximumPenetration_;

   

};


template<typename CAccessor, typename PAccessor>
inline void InitContactsForHCSITS::operator()(size_t c, CAccessor &ca, PAccessor &pa) {

      static_assert(std::is_base_of<data::IContactAccessor, CAccessor>::value, "please provide a valid contact accessor");
      static_assert(std::is_base_of<data::IAccessor, PAccessor>::value, "please provide a valid particle accessor");

   size_t bId1 = ca.getId1(c);
   size_t bId2 = ca.getId2(c);
   ca.setR1(c, ca.getPosition(c) - pa.getPosition(bId1));
   ca.setR2(c, ca.getPosition(c) - pa.getPosition(bId2));


   // construct vector perpendicular to the normal (cross product with cardinal basis vector where the 1 component is where the other vector has its component of smallest magnitude)
   if (std::fabs(ca.getNormalRef(c)[0]) < std::fabs(ca.getNormalRef(c)[1])) {
      if (std::fabs(ca.getNormalRef(c)[0]) < std::fabs(ca.getNormalRef(c)[2]))
         ca.setT(c, Vec3(0, ca.getNormalRef(c)[2], -ca.getNormalRef(c)[1]));   // = n x (1, 0, 0)^T
      else
         ca.setT(c, Vec3(ca.getNormalRef(c)[1], -ca.getNormalRef(c)[0], 0));   // = n x (0, 0, 1)^T
   } else {
      if (std::fabs(ca.getNormalRef(c)[1]) < std::fabs(ca.getNormalRef(c)[2]))
         ca.setT(c, Vec3(-ca.getNormalRef(c)[2], 0, ca.getNormalRef(c)[0]));   // = n x (0, 1, 0)^T
      else
         ca.setT(c, Vec3(ca.getNormalRef(c)[1], -ca.getNormalRef(c)[0], 0));   // = n x (0, 0, 1)^T
   }
   normalize(ca.getTRef(c));
   ca.setO(c, ca.getNormal(c) % ca.getT(c));

   Mat3 contactframe(ca.getNormal(c), ca.getT(c), ca.getO(c));

   // If the distance is negative then penetration is present. This is an error and should be corrected.
   // Correcting the whole error is not recommended since due to the linearization the errors cannot
   // completely fixed anyway and the error reduction will introduce artificial restitution.
   // However, if the distance is positive then it is not about error correction but the distance that
   // can still be overcome without penetration and thus no error correction parameter should be applied.
   if (ca.getDistance(c) < real_t(0.0)) {
      setMaximumPenetration(std::max(getMaximumPenetration(), -ca.getDistance(c)));
      ca.getDistanceRef(c) *= erp_;
   }

   ca.getMuRef(c) = getFriction(pa.getType(bId1), pa.getType(bId2));

   Mat3 diag = -(
           math::skewSymCrossProduct(ca.getR1(c),
157
                                     math::skewSymCrossProduct(pa.getInvInertiaBF(bId1), ca.getR1(c)))
158
           + math::skewSymCrossProduct(ca.getR2(c),
159
                                       math::skewSymCrossProduct(pa.getInvInertiaBF(bId2), ca.getR2(c))));
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
   diag[0] += pa.getInvMass(bId1) + pa.getInvMass(bId2);
   diag[4] += pa.getInvMass(bId1) + pa.getInvMass(bId2);
   diag[8] += pa.getInvMass(bId1) + pa.getInvMass(bId2);

   diag = contactframe.getTranspose() * diag * contactframe;

   // Diagonal block is known to be positive-definite and thus an inverse always exists.
   ca.getDiag_ntoRef(c) = diag;
   ca.getDiag_nto_invRef(c) = diag.getInverse();
   ca.getDiag_n_invRef(c) = math::inv(diag[0]);
   ca.getDiag_to_invRef(c) = Mat2(diag[4], diag[5], diag[7], diag[8]).getInverse();
}

}
}
}