ScaleCoupling.h
Go to the documentation of this file.
1 //Copyright (c) 2013-2023, The MercuryDPM Developers Team. All rights reserved.
2 //For the list of developers, see <http://MercuryDPM.org/Team>.
3 //
4 //Redistribution and use in source and binary forms, with or without
5 //modification, are permitted provided that the following conditions are met:
6 // * Redistributions of source code must retain the above copyright
7 // notice, this list of conditions and the following disclaimer.
8 // * Redistributions in binary form must reproduce the above copyright
9 // notice, this list of conditions and the following disclaimer in the
10 // documentation and/or other materials provided with the distribution.
11 // * Neither the name MercuryDPM nor the
12 // names of its contributors may be used to endorse or promote products
13 // derived from this software without specific prior written permission.
14 //
15 //THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
16 //ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
17 //WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
18 //DISCLAIMED. IN NO EVENT SHALL THE MERCURYDPM DEVELOPERS TEAM BE LIABLE FOR ANY
19 //DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
20 //(INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
21 //LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND
22 //ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
23 //(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
24 //SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
25 
26 #ifndef SCALE_COUPLING_H
27 #define SCALE_COUPLING_H
28 
29 #include <utility>
30 
31 #include "Mercury3D.h"
32 #include "../Coupling/BaseCoupling.h"
33 #include "generic.h"
35 #include "Oomph/SolidProblem.h" //
36 #include "ScaleCoupledElement.h"
37 using namespace oomph;
38 
39 template<class M, class O>
40 class ScaleCoupling : public BaseCoupling<M,O> {
41 
43  typedef typename O::ELEMENT ELEMENT;
44 
48  std::function<double(double x,double y,double z)> couplingWeight_;
49 
51  struct CoupledElement {
52  // pointer to the oomph-lib element
54  // particles in this element
55  std::vector<BaseParticle*> coupledParticles;
56  };
57 
59  std::vector<CoupledElement> coupledElements_;
60 
62  struct CoupledParticle {
66  CoupledElement* coupledElement_pt = nullptr;
68  Vector<double> s;
69 
72  if (coupledElement_pt != nullptr) {
73  remove(coupledElement_pt->coupledParticles.begin(), coupledElement_pt->coupledParticles.end(), particle);
74  coupledElement_pt = nullptr;
75  logger(VERBOSE, "Removed particle % from element %", particle->getIndex(), coupledElement_pt);
76  }
77  }
78 
80  void setCoupledElement(CoupledElement* coupledElementNew, BaseParticle* particle) {
81  if (coupledElement_pt != coupledElementNew) {
82  if (coupledElement_pt != nullptr) {
83  remove(coupledElement_pt->coupledParticles.begin(), coupledElement_pt->coupledParticles.end(), particle);
84  logger(VERBOSE, "Moved particle % from element % to %", particle->getIndex(), coupledElement_pt, coupledElementNew);
85  } else {
86  logger(VERBOSE,"Added particle % to element %", particle->getIndex(), coupledElementNew);
87  }
88  // add coupledElement to particle
89  coupledElement_pt = coupledElementNew;
90  // add particle to coupledElement
91  coupledElement_pt->coupledParticles.push_back(particle);
92  }
93  }
94  };
95 
97  std::vector<CoupledParticle> coupledParticles_;
98 
100  double penalty_ = constants::NaN;
101 
102 public:
103 
105  const std::vector<CoupledElement>& getCoupledElements() {
106  return coupledElements_;
107  }
108 
110  const std::vector<CoupledParticle>& getCoupledParticles() {
111  return coupledParticles_;
112  }
113 
115  void setPenalty(double penalty) {
116  penalty_ = penalty;
117  }
118 
120  void setCouplingWeight(const std::function<double(double x,double y,double z)>& couplingWeight) {
121  couplingWeight_ = couplingWeight;
122  }
123 
126  {
127  // compute nStep
128  unsigned nStep = O::getOomphTimeStep()/M::getTimeStep();
129  if (nStep==0) {
130  // if oomph has a smaller time step than Mercury
131  logger(INFO, "Set nStep %, change mercuryTimeStep from % to %",
132  nStep, M::getTimeStep(), O::getOomphTimeStep());
133  nStep = 1;
134  M::setTimeStep(O::getOomphTimeStep());
135  } else {
136  logger(INFO, "Set nStep %, change oomphTimeStep from % to %",
137  nStep, O::getOomphTimeStep(), nStep * M::getTimeStep());
138  O::setOomphTimeStep(nStep * M::getTimeStep());
139  }
140  // call solve routine
141  solveScaleCoupling(nStep);
142  }
143 
144 private:
145 
147  void solveScaleCoupling(unsigned nStep)
148  {
149  // Check main parameters are set
150  logger.assert_always(O::getOomphTimeStep()>0, "Oomph time step not initialised");
151  logger.assert_always(M::getTimeStep()>0, "Mercury time step not initialised");
152  logger.assert_always(std::isfinite(penalty_), "Penalty parameter not set");
153  logger.assert_always((bool)couplingWeight_, "Coupling weight is not set");
154 
155  // Initialise the coupled element vector
156  initialiseCoupledElements();
157 
158  // This is the main loop advancing over time
159  unsigned nDone = 0; //< last written file number
160  while (M::getTime() < M::getTimeMax())
161  {
162  O::actionsBeforeOomphTimeStep();
163  // solve the coupled problem for one time step
164  computeOneTimeStepScaleCoupling(nStep);
165  // write outputs of the oomphProb; this is slaved to the vtk output of Mercury, i.e. an oomph-lib output get written everytime a Mercury vtk file gets written
166  if (M::getParticlesWriteVTK() && M::getVtkWriter()->getFileCounter() > nDone)
167  {
168  O::writeToVTK();
169  nDone = M::getVtkWriter()->getFileCounter();
170  }
171  }
172 
173  // Write final output files
174  M::finaliseSolve();
175  }
176 
179  // loop through all el_pts, check if any nodes have a fractional coupling weight.
180  // If yes, add the el_pt to the coupledElements_
181  auto nelement = O::mesh_pt()->nelement();
182  for (int e=0; e < nelement; ++e) {
183  bool coupledElement = false;
184  auto* el_pt = dynamic_cast<ELEMENT*>(O::mesh_pt()->finite_element_pt(e));
185  auto nnode = el_pt->nnode();
186  Vector<double> nodalCouplingWeights(nnode);
187  for (int n=0; n < nnode; ++n) {
188  auto* n_pt = dynamic_cast<SolidNode*>(el_pt->node_pt(n));
189  nodalCouplingWeights[n] = couplingWeight_(n_pt->xi(0), n_pt->xi(1), n_pt->xi(2));
190  if (nodalCouplingWeights[n]>0) { // check for nodalCouplingWeights[n]<1 ?
191  coupledElement = true;
192  logger(VERBOSE,"Coupled element at % % %, weight %",n_pt->xi(0), n_pt->xi(1), n_pt->xi(2), nodalCouplingWeights[n]);
193  }
194  }
195  if (coupledElement) {
196  coupledElements_.push_back({el_pt});
197  el_pt->set_coupling_weight(nodalCouplingWeights);
198  }
199  }
200  logger(INFO,"% of % el_pts are coupled", coupledElements_.size(), nelement);
201  }
202 
204  void computeOneTimeStepScaleCoupling(unsigned nStep) {
205  updateCoupledElements();
206  computeCouplingForce();
207  this->solveOomph();
208  this->solveMercury(nStep);
209  };
210 
213  {
215  Vector<double> min {constants::inf, constants::inf, constants::inf};
216  Vector<double> max {-constants::inf, -constants::inf, -constants::inf};
217  // loop over all coupled elements
218  for (auto& coupledElement : coupledElements_)
219  {
220  ELEMENT* el_pt = coupledElement.el_pt;
221  const unsigned nnode = el_pt->nnode();
222  const unsigned dim = el_pt->dim(); // /todo I need to use el_pt->ndim=3, not el_pt->ndim=0, why?
223  // Find min and max
224  for (int n=0; n<nnode; ++n) {
225  auto* n_pt = dynamic_cast<SolidNode *>(el_pt->node_pt(n));
226  for (int d = 0; d < dim; ++d) {
227  min[d] = std::min(min[d], n_pt->x(d));
228  max[d] = std::max(max[d], n_pt->x(d));
229  }
230  }
231  }
232  logger(VERBOSE,"Coupling zone: % < x < %, % < y < %, % < z < % ",min[0],max[0],min[1],max[1],min[2],max[2]);
233 
235  // resize the coupledParticles vector
236  coupledParticles_.resize(M::particleHandler.getSize());
237  // now set the coupledElement pointer for each coupledParticles, and the particles vector for each coupledElement
238  for (int p=0; p<coupledParticles_.size(); ++p) {
239  CoupledParticle& coupledParticle = coupledParticles_[p];
240  Vec3D pos = M::particleHandler.getObject(p)->getPosition();
241  // this is the pointer we need to set
242  if (pos.X>=min[0] && pos.X<=max[0] && pos.Y>=min[1] && pos.Y<=max[1] && pos.Z>=min[2] && pos.Z<=max[2]) {
243  logger(VERBOSE,"Particle % of % is in coupling zone",p, M::particleHandler.getSize());
244  // if near the coupling zone
245  Vector<double> x {pos.X, pos.Y, pos.Z};
246  Vector<double>& s = coupledParticle.s;
247  GeomObject* geom_obj_pt;
248  // if el_pt is already set and still valid, keep
249  if (coupledParticle.coupledElement_pt != nullptr) {
250  coupledParticle.coupledElement_pt->el_pt->locate_zeta(x, geom_obj_pt, s);
251  if (geom_obj_pt!=nullptr) {
252  logger(VERBOSE,"Particle % remains coupled with element %",p, coupledParticle.coupledElement_pt->el_pt);
253  continue;
254  }
255  }
256  // otherwise, search for the right element
257  for (auto& coupledElement : coupledElements_) {
258  coupledElement.el_pt->locate_zeta(x,geom_obj_pt,s);
259  if (geom_obj_pt!=nullptr) {
260  logger(VERBOSE,"Coupled particle % to element % at % % %",p, &coupledElement, s[0], s[1], s[2]);
261  //\todo If a containing element is found, we stop; note though, particles on the boundary could belong to multiple elements; also not clear to me whether force should be applied to global basis function or local
262  coupledParticle.setCoupledElement(&coupledElement, M::particleHandler.getObject(p));
263  break;
264  }
265  }
266  } else {
267  // if away from the coupling, set coupledElement to null
268  if (coupledParticle.coupledElement_pt) logger(VERBOSE, "Removed particle % from element %", p, coupledParticle.coupledElement_pt);
269  coupledParticle.removeCoupledElement(M::particleHandler.getObject(p));
270  }
271  }
272  }
273 
275  Vector<Vector<double>> computeProjectionMatrix(const CoupledElement& coupledElement)
276  {
277  Vector<Vector<double>> projectionMatrix;
278  // get constants
279  const unsigned nParticles = coupledElement.coupledParticles.size();
280  const unsigned nnode = coupledElement.el_pt->nnode();
281  const unsigned dim = coupledElement.el_pt->dim();
282 
283  // prepare projection matrix
284  projectionMatrix.resize(nnode,Vector<double>(nParticles,0.0));
285 
287  // storage for position s and shape psi
288  Vector<double> pos(dim,0.0);
289  Shape psi(nnode);
290  // loop over shape functions at the nodes
291  // loop over shape functions at the particles
292  for (unsigned p = 0; p < nParticles; p++) {
293  BaseParticle* particle = coupledElement.coupledParticles[p];
294  // get the shape function evaluated at the center of particle p
295  coupledElement.el_pt->shape(coupledParticles_[particle->getIndex()].s, psi);
296  for (unsigned n=0; n < nnode; n++) {
297  projectionMatrix[n][p] = psi(n) * particle->getVolume();
298  }
299  }
300  // normalize each row of the projection matrix to sum to one
301  for (unsigned n=0; n < nnode; n++) {
302  double sum = 0;
303  for (unsigned p = 0; p<nParticles; p++) {
304  sum += projectionMatrix[n][p];
305  }
306  if (sum != 0) {
307  for (unsigned p = 0; p < nParticles; p++) {
308  projectionMatrix[n][p] /= sum;
309  }
310  }
311  }
312  return projectionMatrix;
313  }
314 
317  {
318  // first loop over the coupled bulk el_pts
319  for (const auto& coupledElement : coupledElements_)
320  {
321  // construct projection matrix, defined as projectionMatrix_{n,p} = N_{n,p}*V_p / sum_p(N_{n,p}*V_p)
322  // (shape function N, volume V, node n, particle p)
323  Vector<Vector<double>> projectionMatrix = computeProjectionMatrix(coupledElement);
324  // get pointer to the bulk el_pt
325  ELEMENT* el_pt = coupledElement.el_pt;
326  const unsigned nParticles = coupledElement.coupledParticles.size();
327  // get number of nodes in the bulk el_pt
328  const unsigned nnode = coupledElement.el_pt->nnode();
329  // prepare vector of nodal velocities projected from particle centers
330  const unsigned dim = el_pt->dim();
331  // get coupling force at nodal positions, based on penalizing the difference in velocity
332  computeNodalCouplingForces(coupledElement, projectionMatrix);
333  // project coupling forces from nodes to particles (note projectionMatrix^{-1} = projectionMatrix^T)
334  for (unsigned p=0; p < nParticles; p++) {
335  BaseParticle* particle = coupledElement.coupledParticles[p];
336  CoupledParticle& coupledParticle = coupledParticles_[particle->getIndex()];
337  for (unsigned d=0; d < dim; d++) {
338  double forceComponent = 0;
339  for (unsigned n=0; n < nnode; n++) {
340  // note projectionMatrix^{-1} = projectionMatrix^T
341  forceComponent += projectionMatrix[n][p] * el_pt->get_coupling_residual(n, d);
342  }
343  coupledParticle.couplingForce.setComponent(d,forceComponent);
344  }
345  logger(VERBOSE,"Coupling force % on particle %", coupledParticle.couplingForce,p);
346  }
347  }
348  }
349 
351  void computeNodalCouplingForces(const CoupledElement& coupledElement, const Vector<Vector<double>>& projectionMatrix)
352  {
353  // return if no particle in cell
354  if (coupledElement.coupledParticles.empty()) return;
355 
356  // a few shortcut variables
357  const ELEMENT* el_pt = coupledElement.el_pt;
358  const unsigned nnode = el_pt->nnode();
359  const unsigned dim = el_pt->dim();
360  const unsigned nParticles = coupledElement.coupledParticles.size();
361 
363  Vector<Vector<double>> nodalVelocityDifference(nnode,Vector<double>(dim,0.0));
364  // get projected nodal velocity field from particle velocities
365  for (unsigned n=0; n < nnode; n++) {
366  for (unsigned d=0; d < dim; d++) {
367  double nodalVelocityDPM = 0;
368  for (unsigned p=0; p < nParticles; p++) {
369  nodalVelocityDPM += projectionMatrix[n][p] * coupledElement.coupledParticles[p]->getVelocity().getComponent(d);
370  }
371  double nodalVelocityFEM = el_pt->dnodal_position_gen_dt(n,0,d);
372  nodalVelocityDifference[n][d] = nodalVelocityDPM-nodalVelocityFEM;
373  }
374  }
375 
377  Vector<Vector<double> > nodalCouplingForces(nnode, Vector<double>(dim, 0.0));
378 
379  // Set up memory for the shape functions
380  Shape psi(nnode);
381  DShape dpsidxi(nnode,dim);
382  const unsigned n_in = el_pt->integral_pt()->nweight();
383  // Loop over the integration points
384  for (unsigned in=0; in<n_in; in++)
385  {
386  // Get the integral weight
387  double w = el_pt->integral_pt()->weight(in);
388  // Call the derivatives of the shape functions (and get Jacobian)
389  double J = el_pt->dshape_lagrangian_at_knot(in,psi,dpsidxi);
390  // Loop over the test functions, nodes of the element
391  for (unsigned n=0; n < nnode; n++)
392  {
393  // Loop over the nodes of the element again
394  for (unsigned nn=0; nn < nnode; nn++)
395  {
396  double volume = penalty_ * psi(nn) * psi(n) * w * J;
397 
398  for (unsigned d=0; d < dim; d++)
399  {
400  //nodalVelocityDifference[n][d] = 1/penalty_/time_pt()->dt(0);
401  double displacement = nodalVelocityDifference[n][d] * O::time_pt()->dt(0);
402  nodalCouplingForces[n][d] -= volume * displacement;
403  //logger(INFO, "w % J % psi % % displacement % nodalCouplingForces % dt %", w, J, psi(nn), psi(n), displacement, nodalCouplingForces[n][d], time_pt()->dt(0));
404  }
405  }
406  }
407  }
408 
410  coupledElement.el_pt->set_coupling_residual(nodalCouplingForces);
411 
412  logger(VERBOSE,"VelocityDifference % % %, coupling force % % % on element %, node 0, particles %",
413  nodalVelocityDifference[0][0], nodalVelocityDifference[0][1], nodalVelocityDifference[0][2],
414  nodalCouplingForces[0][0], nodalCouplingForces[0][1], nodalCouplingForces[0][2],
415  el_pt, coupledElement.coupledParticles.size());
416  }
417 
420  {
421  if (!p->isFixed())
422  {
423  // Applying the force due to gravity (F = m.g)
424  p->addForce(M::getGravity() * p->getMass());
425  // add particle coupling force 1/w_i*f_i^c to the total force
426  const CoupledParticle& cp = coupledParticles_[p->getIndex()];
427  double couplingWeight = couplingWeight_(p->getPosition().X,p->getPosition().Y,p->getPosition().Z);
428  p->addForce(cp.couplingForce / couplingWeight);
429  }
430  }
431 };
432 
433 #endif //SCALE_COUPLING_H
const unsigned n
Definition: CG3DPackingUnitTest.cpp:32
LL< Log::VERBOSE > VERBOSE
Verbose information.
Definition: Logger.cc:57
LL< Log::INFO > INFO
Info log level.
Definition: Logger.cc:55
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
Definition: BaseCoupling.h:48
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:218
void addForce(const Vec3D &addForce)
Adds an amount to the force on this BaseInteractable.
Definition: BaseInteractable.cc:116
unsigned int getIndex() const
Returns the index of the object in the handler.
Definition: BaseObject.h:118
Definition: BaseParticle.h:54
bool isFixed() const override
Is fixed Particle function. It returns whether a Particle is fixed or not, by checking its inverse Ma...
Definition: BaseParticle.h:93
virtual Mdouble getVolume() const
Get Particle volume function, which required a reference to the Species vector. It returns the volume...
Definition: BaseParticle.cc:143
Mdouble getMass() const
Returns the particle's mass.
Definition: BaseParticle.h:322
Definition: ScaleCoupling.h:40
void solveScaleCoupling(unsigned nStep)
Solve scale-coupled problem.
Definition: ScaleCoupling.h:147
void setCouplingWeight(const std::function< double(double x, double y, double z)> &couplingWeight)
Sets weight function (determines which nodes/el_pts are in coupling zone)
Definition: ScaleCoupling.h:120
void computeCouplingForce()
Computes coupling force for each element and particle.
Definition: ScaleCoupling.h:316
void updateCoupledElements()
Updates the coupledElements and coupledParticles vectors: finds which particles interact with which e...
Definition: ScaleCoupling.h:212
O::ELEMENT ELEMENT
type of el_pt in O (e.g. RefineableQDPVDElement<3, 2>)
Definition: ScaleCoupling.h:43
std::function< double(double x, double y, double z)> couplingWeight_
Definition: ScaleCoupling.h:48
std::vector< CoupledElement > coupledElements_
Vector of all coupled elements.
Definition: ScaleCoupling.h:59
const std::vector< CoupledParticle > & getCoupledParticles()
get coupled particle
Definition: ScaleCoupling.h:110
void computeOneTimeStepScaleCoupling(unsigned nStep)
What to do in each time step.
Definition: ScaleCoupling.h:204
const std::vector< CoupledElement > & getCoupledElements()
get coupled element
Definition: ScaleCoupling.h:105
void solveScaleCoupling()
Computes nStep, the ratio of FEM and DEM time steps, then calls solveSurfaceCoupling(nStep)
Definition: ScaleCoupling.h:125
std::vector< CoupledParticle > coupledParticles_
The i-th element of this vector describes the coupling properties of the i-th DPM particle.
Definition: ScaleCoupling.h:97
void computeExternalForces(BaseParticle *p) override
Applies coupling force to MercuryDPM in each time step.
Definition: ScaleCoupling.h:419
void computeNodalCouplingForces(const CoupledElement &coupledElement, const Vector< Vector< double >> &projectionMatrix)
Computes coupling force for each node.
Definition: ScaleCoupling.h:351
void setPenalty(double penalty)
Sets penalty parameter.
Definition: ScaleCoupling.h:115
void initialiseCoupledElements()
initialises the coupledElements vector
Definition: ScaleCoupling.h:178
Vector< Vector< double > > computeProjectionMatrix(const CoupledElement &coupledElement)
Computes projectionMatrix.
Definition: ScaleCoupling.h:275
Definition: Vector.h:51
Mdouble Y
Definition: Vector.h:66
Mdouble Z
Definition: Vector.h:66
Mdouble X
the vector components
Definition: Vector.h:66
void setComponent(int index, double val)
Sets the requested component of this Vec3D to the requested value.
Definition: Vector.cc:217
const Mdouble NaN
Definition: GeneralDefine.h:43
const Mdouble inf
Definition: GeneralDefine.h:44
Definition: AnisotropicHookean.h:31
Stores properties of a coupled element: pointer to the element and list of particles coupled to it.
Definition: ScaleCoupling.h:51
std::vector< BaseParticle * > coupledParticles
Definition: ScaleCoupling.h:55
ELEMENT * el_pt
Definition: ScaleCoupling.h:53
For all particles, stores coupling properties: coupling force, pointer to coupled element and locatio...
Definition: ScaleCoupling.h:62
CoupledElement * coupledElement_pt
Element in which this particle resides (null if particle is not coupled)
Definition: ScaleCoupling.h:66
Vec3D couplingForce
Coupling force.
Definition: ScaleCoupling.h:64
void setCoupledElement(CoupledElement *coupledElementNew, BaseParticle *particle)
Remove particle from coupledParticle.coupledElement.particles and set coupledParticle....
Definition: ScaleCoupling.h:80
void removeCoupledElement(BaseParticle *particle)
Removes particle from coupledParticle.coupledElement.particles and sets coupledParticle....
Definition: ScaleCoupling.h:71
Vector< double > s
Location of particle in coupled element, returned by locate zeta.
Definition: ScaleCoupling.h:68