Alamo
Fracture.H
Go to the documentation of this file.
1 #ifndef INTEGRATOR_FRACTURE_H
2 #define INTEGRATOR_FRACTURE_H
3 
4 #include <iostream>
5 #include <fstream>
6 #include <iomanip>
7 
8 #include "AMReX.H"
9 #include "AMReX_ParallelDescriptor.H"
10 #include "AMReX_ParmParse.H"
11 
12 #include "Integrator/Integrator.H"
13 
14 #include "BC/BC.H"
15 #include "BC/Constant.H"
16 
17 #include "IC/IC.H"
18 #include "IC/Ellipsoid.H"
19 #include "IC/Notch.H"
20 #include "IC/Laminate.H"
21 #include "IC/PerturbedInterface.H"
22 
23 #include "Operator/Operator.H"
24 #include "Operator/Elastic.H"
25 #include "Solver/Nonlocal/Linear.H"
26 #include "Solver/Nonlocal/Newton.H"
27 
28 #include "Model/Solid/Solid.H"
30 //#include "Model/Solid/Linear/IsotropicDegradable.H"
31 //#include "Model/Solid/Affine/J2PlasticDegradable.H"
32 //#include "Model/Solid/Affine/CrystalPlasticDegradable.H"
33 
37 
38 #include "Numeric/Stencil.H"
39 #include <eigen3/Eigen/Dense>
40 
41 namespace Integrator
42 {
44 //using brittle_fracture_model_type_test = Model::Solid::Linear::IsotropicDegradable;
45 
46 class Fracture : public Integrator
47 {
48 
49 public:
51  {
52  IO::ParmParse pp_crack("crack");
53  pp_crack.query("modulus_scaling_max",crack.scaleModulusMax);
54  pp_crack.query("refinement_threshold",crack.refinement_threshold);
55 
56 
57  pp_crack.queryclass("constant",crack.cracktype);
58 
59  IO::ParmParse pp_crack_df("crack.df");
60  pp_crack_df.query("mult_Gc", crack.mult_df_Gc);
61  pp_crack_df.query("mult_Lap", crack.mult_df_lap);
62  pp_crack_df.query("beta", crack.beta);
63  pp_crack_df.query("tol_rel",crack.driving_force_tolerance_rel);
64  pp_crack_df.query("tol_abs",crack.driving_force_tolerance_abs);
65 
66  // IC for crack field
67  IO::ParmParse pp("crack.ic");
68  pp_queryarr("type", crack.ic_type); // Type of crack {notch,ellipsoid}
69 
70  if (crack.ic_type.size() == 0) crack.is_ic = false;
71  else
72  {
73  crack.is_ic = true;
74  for (int i = 0; i < crack.ic_type.size(); i++)
75  {
76  if(crack.ic_type[i] == "notch")
77  {
78  IC::Notch *tmpic = new IC::Notch(geom);
79  pp_queryclass("notch",*tmpic);
80  crack.ic.push_back(tmpic);
81  }
82  else if (crack.ic_type[i] == "ellipsoid")
83  {
84  IC::Ellipsoid *tmpic = new IC::Ellipsoid(geom);
85  pp_queryclass("ellipsoid",*tmpic);
86  crack.ic.push_back(tmpic);
87  }
88  else
89  Util::Abort(INFO, "This IC hasn't been implemented yet");
90  }
91  }
92 
93  RegisterNodalFab(crack.c_mf, 1, number_of_ghost_nodes, "crack", true);
94  RegisterNodalFab(crack.c_old_mf, 1, number_of_ghost_nodes, "crack_old", false);
95  RegisterNodalFab(crack.driving_force_mf, 5, number_of_ghost_nodes, "driving_force", false);
96 
97  RegisterIntegratedVariable(&(crack.driving_force_norm),"driving_force_norm");
98  RegisterIntegratedVariable(&(crack.int_crack),"int_c");
99  RegisterIntegratedVariable(&(elastic.int_energy),"int_W");
100 
101  IO::ParmParse pp_material("material");
102  pp_material.query("refinement_threshold",material.refinement_threshold);
103 
104  IO::ParmParse pp_material_ic("material.ic");
105  pp_material_ic.query("type",material.ic_type);
106  if (material.ic_type == "laminate")
107  {
108  IC::Laminate *tmpic = new IC::Laminate(geom);
109  pp_material_ic.queryclass("laminate",*tmpic);
110 
111  IO::ParmParse _temp("material.ic.laminate");
112  _temp.query("number_of_inclusions",number_of_materials);
114 
115  material.ic = tmpic;
116  material.is_ic = true;
117  }
118  else if (material.ic_type == "ellipse")
119  {
120  IC::Ellipse *tmpic = new IC::Ellipse(geom);
121  pp_material_ic.queryclass("ellipse",*tmpic);
123 
124  material.ic = tmpic;
125  material.is_ic = true;
126  }
127  else if (material.ic_type == "perturbed_interface")
128  {
130  pp_material_ic.queryclass("perturbed_interface", *tmpic);
132 
133  material.ic = tmpic;
134  material.is_ic = true;
135  }
136  else
137  material.is_ic = false;
138 
139  for (int p = 0; p < number_of_materials; p++)
140  {
141  std::string pp_string = "material.model"+std::to_string(p+1);
142  Util::Message(INFO, pp_string);
143  IO::ParmParse pp_model(pp_string);
144 
145  if(pp_model.contains("type"))
146  {
147  std::string model_type_str;
148  pp_model.query("type", model_type_str);
149  if(model_type_str == "isotropic")
150  {
152  pp_model.queryclass("isotropic",*tmpmodel);
153  material.brittlemodeltype.push_back(*tmpmodel);
154  }
155  }
156  else
157  {
158  if (p == 0) Util::Abort(INFO, "Need material information for at least one material");
159  material.brittlemodeltype.push_back(material.brittlemodeltype[p-1]);
160  }
161  }
162 
163  IO::ParmParse pp_load("loading");
164  if (pp_load.countval("body_force")) pp_load.queryarr("body_force",loading.body_force);
165  pp_load.query("val", loading.val);
166 
167  IO::ParmParse pp_elastic("elastic");
168  pp_elastic.query("df_mult", elastic.df_mult);
169  pp_elastic.queryclass("bc",elastic.brittlebc);
170  pp_elastic.query("int",elastic.interval);
171  pp_elastic.query("omega",elastic.omega);
172 
173 
174  nlevels = maxLevel() + 1;
175  RegisterNodalFab(elastic.disp_mf, AMREX_SPACEDIM, number_of_ghost_nodes, "disp", true);
176  RegisterNodalFab(elastic.rhs_mf, AMREX_SPACEDIM, number_of_ghost_nodes, "rhs", false);
177  RegisterNodalFab(elastic.residual_mf, AMREX_SPACEDIM, number_of_ghost_nodes, "res", false);
178  RegisterNodalFab(elastic.strain_mf, AMREX_SPACEDIM*AMREX_SPACEDIM, number_of_ghost_nodes, "strain", false);
179  RegisterNodalFab(elastic.stress_mf, AMREX_SPACEDIM*AMREX_SPACEDIM, number_of_ghost_nodes, "stress", true);
180  RegisterNodalFab(elastic.energy_mf, 1, number_of_ghost_nodes, "energy", false);
181  RegisterNodalFab(elastic.energy_pristine_mf, 1, number_of_ghost_nodes, "energy_pristine", false);
182  RegisterNodalFab(elastic.energy_pristine_old_mf, 1, number_of_ghost_nodes, "energy_pristine_old", false);
183  RegisterNodalFab(material.material_mf, number_of_materials, number_of_ghost_nodes, "materials", true);
185 
186  }
187 
188 protected:
189  void Initialize(int ilev) override
190  {
191  elastic.disp_mf[ilev]->setVal(0.0);
192  elastic.strain_mf[ilev]->setVal(0.0);
193  elastic.stress_mf[ilev]->setVal(0.0);
194  elastic.rhs_mf[ilev]->setVal(0.0);
195  elastic.energy_mf[ilev]->setVal(0.0);
196  elastic.residual_mf[ilev]->setVal(0.0);
197  elastic.energy_pristine_mf[ilev] -> setVal(0.);
198  elastic.energy_pristine_old_mf[ilev] -> setVal(0.);
199 
200  if(crack.is_ic)
201  {
202  crack.c_mf[ilev]->setVal(1.0);
203  crack.c_old_mf[ilev]->setVal(1.0);
204 
205  for (int i = 0; i < crack.ic.size(); i++)
206  {
207  crack.ic[i]->Add(ilev,crack.c_mf);
208  crack.ic[i]->Add(ilev,crack.c_old_mf);
209  }
210  }
211  else
212  {
213  crack.c_mf[ilev]->setVal(1.0);
214  crack.c_old_mf[ilev]->setVal(1.0);
215  }
216 
217  if(material.is_ic) material.ic->Initialize(ilev,material.material_mf);
218  else material.material_mf[ilev]->setVal(1.0);
219  }
220 
221  void TimeStepBegin(Set::Scalar /*time*/, int /*iter*/) override
222  {
223  Util::Message(INFO,crack.driving_force_norm," ",crack.driving_force_reference," ",crack.driving_force_tolerance_rel);
224  if (crack.driving_force_norm / crack.driving_force_reference < crack.driving_force_tolerance_rel)
225  elastic.do_solve_now = true;
226  if (crack.driving_force_norm < crack.driving_force_tolerance_abs)
227  elastic.do_solve_now = true;
228 
229  if (!elastic.do_solve_now) return;
230  // if(iter%elastic.interval != 0) return;
231 
233  Util::Abort(INFO,"DegradeModulus is no longer implemented");
234  // for (int p = 0; p < number_of_materials; p++)
235  // material.brittlemodeltype[p].DegradeModulus(0.0);
236 
237  for (int ilev = 0; ilev < nlevels; ++ilev)
238  {
239  //material.model_mf[ilev]->setVal(material.brittlemodeltype[0]);
240  Util::RealFillBoundary(*material.material_mf[ilev],geom[ilev]);
241  std::swap(elastic.energy_pristine_old_mf[ilev], elastic.energy_pristine_mf[ilev]);
242 
243  for (MFIter mfi(*material.material_mf[ilev], false); mfi.isValid(); ++mfi)
244  {
245  amrex::Box bx = mfi.grownnodaltilebox();
246  amrex::Array4<brittle_fracture_model_type_test> const &model = material.model_mf[ilev]->array(mfi);
247  amrex::Array4<const Set::Scalar> const &mat = material.material_mf[ilev]->array(mfi);
248 
249  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
250  brittle_fracture_model_type_test modelsum(0.0,0.0);
251  for (int p = 0; p < number_of_materials; p++)
252  {
253  modelsum += mat(i,j,k,p)*material.brittlemodeltype[p];
254  }
255  model(i,j,k) = /*(1./matsum) * */modelsum ;
256  });
257  }
258  Util::RealFillBoundary(*material.model_mf[ilev],geom[ilev]);
259  }
260 
261  // Degrading the modulus
262  for (int ilev = 0; ilev < nlevels; ++ilev)
263  {
264  Util::RealFillBoundary(*crack.c_mf[ilev],geom[ilev]);
265  Util::RealFillBoundary(*material.material_mf[ilev],geom[ilev]);
266 
267  for (amrex::MFIter mfi(*elastic.disp_mf[ilev],true); mfi.isValid(); ++mfi)
268  {
269  amrex::Box box = mfi.grownnodaltilebox();
270  amrex::Array4<const Set::Scalar> const& c = (*crack.c_mf[ilev]).array(mfi);
271  Util::Abort(INFO,"This is depricated and needs to be upgraded!");
272  //amrex::Array4<brittle_fracture_model_type_test> model = (material.model_mf)[ilev]->array(mfi);
273 
274  amrex::ParallelFor (box,[=] AMREX_GPU_DEVICE(int i, int j, int k){
275  Set::Scalar _temp = 0.0;
276  // _temp = crack.cracktype->g_phi(c(i,j,k,0),0);
277  _temp = crack.cracktype.g_phi(c(i,j,k,0),0) + crack.scaleModulusMax;
278 
279  if (std::isnan(_temp)) Util::Abort(INFO);
280  if(_temp < 0.0) _temp = 0.;
281  if(_temp > 1.0) _temp = 1.0;
282 
283  // _temp = std::max(crack.scaleModulusMax, _temp);
284  // _temp = std::max(crack.scaleModulusMax, crack.scaleModulusMax + _temp * (1. - crack.scaleModulusMax));
285 
286  Util::Abort(INFO,"DegradeModulus is no longer available");
287  //model(i,j,k).DegradeModulus(1-_temp);
288  });
289  }
290  Util::RealFillBoundary(*material.model_mf[ilev],geom[ilev]);
291  }
292 
293  //
294  // Update the body force to be consistent with grid size
295  //
296  for (int ilev = 0; ilev < nlevels; ++ilev)
297  {
298  const Real* DX = geom[ilev].CellSize();
299  Set::Scalar volume = AMREX_D_TERM(DX[0],*DX[1],*DX[2]);
300 
301  AMREX_D_TERM(elastic.rhs_mf[ilev]->setVal(loading.body_force(0)*volume,0,1);,
302  elastic.rhs_mf[ilev]->setVal(loading.body_force(1)*volume,1,1);,
303  elastic.rhs_mf[ilev]->setVal(loading.body_force(2)*volume,2,1););
304  }
305 
306  elastic.brittlebc.Init(elastic.rhs_mf,geom);
307 
308  //
309  // Prepare the operator for solve
310  //
312  LPInfo info;
313  for (int ilev = 0; ilev < nlevels; ilev++) if (elastic.disp_mf[ilev]->contains_nan()) Util::Warning(INFO);
314 
315  op_b.define(geom, grids, dmap, info);
316  op_b.SetBC(&elastic.brittlebc);
317  op_b.SetAverageDownCoeffs(true); // <<< Important! Do not change.
318  op_b.SetOmega(elastic.omega);
319  //
320  // Actually do the elastic solve
321  //
322  //for (int ilev = 0; ilev < nlevels; ilev++)
323  //{
324  // elastic.disp_mf[ilev]->setVal(0.0);
325  //}
326 
327  {
328  IO::ParmParse pp("elastic");
330  //Solver::Nonlocal::Linear<brittle_fracture_model_type_test> solver(op_b);
331  pp_queryclass("solver",solver);
332  solver.solve(elastic.disp_mf, elastic.rhs_mf, material.model_mf,1E-8,1E-8);
333  solver.compResidual(elastic.residual_mf,elastic.disp_mf,elastic.rhs_mf,material.model_mf);
334  }
335 
336  for (int ilev = 0; ilev < nlevels; ilev++)
337  {
338  op_b.Strain(ilev,*elastic.strain_mf[ilev],*elastic.disp_mf[ilev]);
339  op_b.Stress(ilev,*elastic.stress_mf[ilev],*elastic.disp_mf[ilev]);
340  op_b.Energy(ilev,*elastic.energy_mf[ilev],*elastic.disp_mf[ilev]);
341 
342  Util::RealFillBoundary(*elastic.strain_mf[ilev],geom[ilev]);
343  Util::RealFillBoundary(*elastic.energy_pristine_old_mf[ilev],geom[ilev]);
344  Util::RealFillBoundary(*material.material_mf[ilev],geom[ilev]);
345  elastic.energy_pristine_mf[ilev]->setVal(0.0);
346 
347  for (amrex::MFIter mfi(*elastic.strain_mf[ilev],true); mfi.isValid(); ++mfi)
348  {
349  const amrex::Box& box = mfi.grownnodaltilebox();
350  amrex::Array4<const Set::Scalar> const& strain = (*elastic.strain_mf[ilev]).array(mfi);
351  amrex::Array4<const Set::Scalar> const& mat = (*material.material_mf[ilev]).array(mfi);
352  amrex::Array4<Set::Scalar> const& energy = (*elastic.energy_pristine_mf[ilev]).array(mfi);
353  amrex::Array4<Set::Scalar> const& energy_old = (*elastic.energy_pristine_old_mf[ilev]).array(mfi);
354 
355  amrex::ParallelFor (box,[=] AMREX_GPU_DEVICE(int i, int j, int k){
356 
357  Set::Matrix eps = Numeric::FieldToMatrix(strain,i,j,k);
358  Eigen::SelfAdjointEigenSolver<Set::Matrix> eigensolver(eps);
359  Set::Vector eValues = eigensolver.eigenvalues();
360  Set::Matrix eVectors = eigensolver.eigenvectors();
361 
362  Set::Matrix epsp = Set::Matrix::Zero();
363  Set::Matrix epsn = Set::Matrix::Zero();
364 
365  for (int n = 0; n < AMREX_SPACEDIM; n++)
366  {
367  if(eValues(n) > 0.0) epsp += eValues(n)*(eVectors.col(n)*eVectors.col(n).transpose());
368  else epsn += eValues(n)*(eVectors.col(n)*eVectors.col(n).transpose());
369  }
370  // brittle_fracture_model_type_test _tmpmodel(0.0,0.0);
371  // for (int p = 0; p < number_of_materials; p++)
372  // _tmpmodel += mat(i,j,k,p)*material.brittlemodeltype[p];
373 
374 
375  for (int p = 0; p < number_of_materials; p++)
376  energy(i,j,k) += mat(i,j,k,p)*material.brittlemodeltype[p].W(epsp);
377 
378  if (energy(i,j,k) < energy_old(i,j,k)) energy(i,j,k) = energy_old(i,j,k);
379  });
380  }
381  Util::RealFillBoundary(*elastic.energy_pristine_mf[ilev],geom[ilev]);
382  }
383 
386  }
387 
388  void Advance(int lev, Set::Scalar /*time*/, Set::Scalar dt) override
389  {
390  std::swap(crack.c_old_mf[lev], crack.c_mf[lev]);
391 
392  //Util::RealFillBoundary(*crack.c_old_mf[lev],geom[lev]);
393  //Util::RealFillBoundary(*crack.c_mf[lev],geom[lev]);
394  //Util::RealFillBoundary(*material.material_mf[lev],geom[lev]);
395  //Util::RealFillBoundary(*elastic.energy_pristine_mf[lev],geom[lev]);
396 
397  const Set::Scalar* DX = geom[lev].CellSize();
398  amrex::Box domain(geom[lev].Domain());
399  domain.convert(amrex::IntVect::TheNodeVector());
400  const amrex::Dim3 lo= amrex::lbound(domain), hi = amrex::ubound(domain);
401 
402  for ( amrex::MFIter mfi(*crack.c_mf[lev],true); mfi.isValid(); ++mfi )
403  {
404  const amrex::Box& bx = mfi.nodaltilebox();
405  amrex::Array4<const Set::Scalar> const& c_old = (*crack.c_old_mf[lev]).array(mfi);
406  amrex::Array4<Set::Scalar> const& df = (*crack.driving_force_mf[lev]).array(mfi);
407  amrex::Array4<Set::Scalar> const& c = (*crack.c_mf[lev]).array(mfi);
408  amrex::Array4<const Set::Scalar> const& energy = (*elastic.energy_pristine_mf[lev]).array(mfi);
409  amrex::Array4<const Set::Scalar> const& mat = (*material.material_mf[lev]).array(mfi);
410 
411  amrex::ParallelFor (bx,[=] AMREX_GPU_DEVICE(int i, int j, int k){
412 
413 #if AMREX_SPACEDIM !=2
414  Util::Abort(INFO, "This doesn't work with 1D or 3D yet");
415 #endif
416  if (i == lo.x && j == lo.y) c(i,j,k,0) = c(i+1,j+1,k,0);
417  else if (i == lo.x && j == hi.y) c(i,j,k,0) = c(i+1,j-1,k,0);
418  else if (i == hi.x && j == lo.y) c(i,j,k,0) = c(i-1,j+1,k,0);
419  else if (i == hi.x && j == hi.y) c(i,j,k,0) = c(i-1,j-1,k,0);
420  else if (i == lo.x) c(i,j,k,0) = c(i+1,j,k,0);
421  else if (j == lo.y) c(i,j,k,0) = c(i,j+1,k,0);
422  else if (i == hi.x) c(i,j,k,0) = c(i-1,j,k,0);
423  else if (j == hi.y) c(i,j,k,0) = c(i,j-1,k,0);
424  else
425  {
426  Set::Scalar rhs = 0.0;
427  //Set::Vector Dc = Numeric::Gradient(c_old, i, j, k, 0, DX);
428 
429  Set::Scalar bilap =
430  Numeric::Stencil<Set::Scalar, 4, 0, 0>::D(c_old,i,j,k,0,DX) +
431  2.0 * Numeric::Stencil<Set::Scalar, 2, 2, 0>::D(c_old,i,j,k,0,DX) +
433 
434  Set::Scalar en_cell = energy(i,j,k,0);
435 
436  //if (std::isnan(en_cell)) Util::Abort(INFO, "Nans detected in en_cell. energy_box(i,j,k,0) = ", energy(i,j,k,0));
437  //if (std::isinf(c_old(i,j,k,0))) Util::Abort(INFO, "Infs detected in c_old");
438  //if (c_old(i,j,k,0) > 1) Util::Abort(INFO, "Very large values of c_old at (",i,",",j,",",k,") = ", c_old(i,j,k,0));
439 
440  df(i,j,k,0) = crack.cracktype.Dg_phi(c_old(i,j,k,0),0.0)*en_cell*elastic.df_mult;
441  rhs += crack.cracktype.Dg_phi(c_old(i,j,k,0),0.0)*en_cell*elastic.df_mult;
442 
443  //if(rhs > 1.e10) Util::Abort(INFO, "Very large values of rhs at (",i,",",j,",",k,") = ", rhs);
444 
445  //Set::Matrix DDc = Numeric::Hessian(c_old, i, j, k, 0, DX);
446  Set::Scalar laplacian = Numeric::Laplacian(c_old,i,j,k,0,DX);//DDc.trace();
447 
448  Set::Scalar Gc = crack.cracktype.Gc(0.0);
449  Set::Scalar Zeta = crack.cracktype.Zeta(0.0);
450  Set::Scalar _temp_product = 1.0;
451  for (int m = 0; m < number_of_materials; m++)
452  _temp_product *= mat(i,j,k,m);
453 
454  if (number_of_materials > 1) Gc *= (1.0 - _temp_product*(1.-crack.mult_df_Gc)/0.25);
455 
456  // Util::Message(INFO, "Gc = ", Gc);
457 
458  df(i,j,k,1) = Gc*crack.cracktype.Dw_phi(c_old(i,j,k,0),0.0)/(4.0*Zeta);
459  rhs += Gc*crack.cracktype.Dw_phi(c_old(i,j,k,0),0.)/(4.0*Zeta);
460 
461  df(i,j,k,2) = 2.0*Zeta*Gc*laplacian*crack.mult_df_lap;
462  rhs -= 2.0*Zeta*Gc*laplacian*crack.mult_df_lap;
463 
464  df(i,j,k,3) = crack.beta*bilap;
465  rhs += crack.beta * bilap;
466 
467  //if (!material.homogeneous)
468  //rhs *= crack.cracktype->g_phi(_temp,0.0);
469 
470  //if(std::isnan(rhs)) Util::Abort(INFO, "Dwphi = ", crack.cracktype->Dw_phi(c_old(i,j,k,0),0.0),". c_old(i,j,k,0) = ",c_old(i,j,k,0));
471  df(i,j,k,4) = std::max(0.,rhs - crack.cracktype.DrivingForceThreshold(c_old(i,j,k,0)));
472  c(i,j,k,0) = c_old(i,j,k,0) - dt*std::max(0.,rhs - crack.cracktype.DrivingForceThreshold(c_old(i,j,k,0)))*crack.cracktype.Mobility(c_old(i,j,k,0));
473 
474  if (c (i,j,k,0) < 0.0) c(i,j,k,0) = 0.0;
475  if (c (i,j,k,0) > 1.0) c(i,j,k,0) = 1.0;
476  }
477  });
478  }
479  Util::RealFillBoundary(*crack.c_mf[lev],geom[lev]);
480  }
481 
482  void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, amrex::Real /*time*/, int /*ngrow*/) override
483  {
484  const amrex::Real *DX = geom[lev].CellSize();
485  const Set::Vector dx(DX);
486  const Set::Scalar dxnorm = dx.lpNorm<2>();
487  amrex::Box domain(geom[lev].Domain());
488  domain.convert(amrex::IntVect::TheNodeVector());
489 
490  for (amrex::MFIter mfi(*crack.c_mf[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
491  {
492  const amrex::Box bx = mfi.tilebox() & domain.grow(-1);
493  amrex::Array4<char> const &tags = a_tags.array(mfi);
494  amrex::Array4<const Set::Scalar> const &c = (*crack.c_mf[lev]).array(mfi);
495  amrex::Array4<const Set::Scalar> const &mat = (*material.material_mf[lev]).array(mfi);
496 
497  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
498  Set::Vector grad = Numeric::Gradient(c, i, j, k, 0, DX);
499  Set::Vector grad2 = Numeric::Gradient(mat, i, j, k, 0, DX);
500  if (dxnorm * grad.lpNorm<2>() >= crack.refinement_threshold || dxnorm * grad2.lpNorm<2>() >= material.refinement_threshold)
501  tags(i, j, k) = amrex::TagBox::SET;
502  });
503  }
504  }
505 
506  void Integrate(int amrlev, Set::Scalar /*time*/, int /*step*/,const amrex::MFIter &mfi, const amrex::Box &box) override
507  {
508  const amrex::Real* DX = geom[amrlev].CellSize();
509  const Set::Scalar DV = AMREX_D_TERM(DX[0],*DX[1],*DX[2]);
510 
511  amrex::Array4<const Set::Scalar> const &df = (*crack.driving_force_mf[amrlev]).array(mfi);
512  amrex::Array4<const Set::Scalar> const &c_new = (*crack.c_mf[amrlev]).array(mfi);
513  amrex::Array4<const Set::Scalar> const &energy = (*elastic.energy_mf[amrlev]).array(mfi);
514 
515  amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
516  {
517  crack.driving_force_norm += df(i,j,k,4) * DV;
518  crack.int_crack += c_new(i,j,k) * DV;
519  elastic.int_energy += energy(i,j,k) * DV;
520  });
521  }
522 
523  void TimeStepComplete(amrex::Real /*time*/,int /*iter*/) override
524  {
525  if (elastic.do_solve_now)
526  crack.driving_force_reference = crack.driving_force_norm;
527  //crack.driving_force_reference = crack.driving_force_norminf;
528 
529  elastic.do_solve_now = false;
530  }
531 
532 private:
533 
534  int number_of_ghost_nodes = 2; ///< Number of ghost nodes
536  int nlevels;
537 
538  struct{
539  Set::Field<Set::Scalar> disp_mf; ///< displacement field
540  Set::Field<Set::Scalar> strain_mf; ///< total strain field (gradient of displacement)
542  Set::Field<Set::Scalar> rhs_mf; ///< rhs fab for elastic solution
543  Set::Field<Set::Scalar> residual_mf; ///< residual field for solver
544  Set::Field<Set::Scalar> energy_mf; ///< total elastic energy
545  Set::Field<Set::Scalar> energy_pristine_mf; ///< energy of the prisitne material as if no crack is present
546  Set::Field<Set::Scalar> energy_pristine_old_mf; ///< energy of the pristine material for previous time step.
547  Set::Scalar int_energy = 0.0; ///< integrated energy over the entire domain.
548 
549  BC::Operator::Elastic::Constant brittlebc; ///< elastic BC if using brittle fracture
550  Set::Scalar df_mult = 1.0; ///< mulitplier for elastic driving force.
551  bool do_solve_now = false;
552  int interval = 0;
554  } elastic;
555 
556 
557 
558  struct{
559  Set::Field<Set::Scalar> c_mf; ///< crack field at current time step
560  Set::Field<Set::Scalar> c_old_mf; ///< crack field at previous time step
561  Set::Field<Set::Scalar> driving_force_mf; ///< crack driving forces.
562  Set::Scalar int_crack = 0.0; ///< integrated crack field over the entire domain.
568 
569  Model::Interface::Crack::Constant cracktype; ///< type of crack. See Crack/Constant or Crack/Sin
570  amrex::Vector<std::string> ic_type; ///< crack IC type. See IC/Notch and IC/Ellipsoid
571  amrex::Vector<IC::IC*> ic; ///< crack IC. See IC/Notch and IC/Ellipsoid
572  bool is_ic = false;
573 
574  Set::Scalar scaleModulusMax = 0.02; ///< material modulus ratio inside crack (default = 0.02).
576 
577  Set::Scalar mult_df_Gc = 1.0; ///< Multiplier for Gc/zeta term
578  Set::Scalar mult_df_lap = 1.0; ///< Multiplier for the laplacian term
580  } crack;
581 
582  struct{
583  amrex::Vector<brittle_fracture_model_type_test> brittlemodeltype;
586  std::string input_material = "isotropic";
588 
590  std::string ic_type;
591  bool is_ic = false;
592  } material;
593 
594  struct{
595  Set::Vector body_force = Set::Vector::Zero();
597  } loading;
598 
599 };
600 }
601 #endif
Integrator::Fracture::Advance
void Advance(int lev, Set::Scalar, Set::Scalar dt) override
Definition: Fracture.H:388
Integrator::Fracture::ic_type
amrex::Vector< std::string > ic_type
crack IC type. See IC/Notch and IC/Ellipsoid
Definition: Fracture.H:570
Util::RealFillBoundary
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T >> &a_mf, const amrex::Geometry &)
Definition: Util.H:307
Integrator::Fracture::crack
struct Integrator::Fracture::@8 crack
Integrator::Fracture::driving_force_norm
Set::Scalar driving_force_norm
Definition: Fracture.H:564
Integrator::Fracture::driving_force_tolerance_abs
Set::Scalar driving_force_tolerance_abs
Definition: Fracture.H:567
IC::IC
Pure abstract IC object from which all other IC objects inherit.
Definition: IC.H:12
IC::Ellipse
Definition: Ellipse.H:21
Operator::Elastic
Definition: Elastic.H:22
Integrator::Fracture::is_ic
bool is_ic
Definition: Fracture.H:572
Integrator::Integrator::box
std::vector< amrex::Box > box
Definition: Integrator.H:367
PerturbedInterface.H
Set::Isotropic
@ Isotropic
Definition: Base.H:197
Integrator::Fracture
Definition: Fracture.H:46
Integrator::Fracture::disp_mf
Set::Field< Set::Scalar > disp_mf
displacement field
Definition: Fracture.H:539
Integrator.H
Notch.H
Set::Field< Set::Scalar >
Definition: Set.H:236
Integrator::Fracture::driving_force_mf
Set::Field< Set::Scalar > driving_force_mf
crack driving forces.
Definition: Fracture.H:561
Integrator::Integrator::integrate_variables_after_advance
bool integrate_variables_after_advance
Definition: Integrator.H:314
Laminate.H
Solver::Nonlocal::Newton
Definition: Newton.H:15
Set::Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition: Base.H:20
Constant.H
Numeric::Gradient
AMREX_FORCE_INLINE Set::Vector Gradient(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const int &m, const Set::Scalar dx[AMREX_SPACEDIM], std::array< StencilType, AMREX_SPACEDIM > stencil=DefaultType)
Definition: Stencil.H:655
Integrator::Fracture::model_mf
Set::Field< brittle_fracture_model_type_test > model_mf
Definition: Fracture.H:585
Elastic.H
Integrator::Fracture::Integrate
void Integrate(int amrlev, Set::Scalar, int, const amrex::MFIter &mfi, const amrex::Box &box) override
Definition: Fracture.H:506
Integrator::Fracture::number_of_ghost_nodes
int number_of_ghost_nodes
Number of ghost nodes.
Definition: Fracture.H:534
Integrator::Fracture::Initialize
void Initialize(int ilev) override
Definition: Fracture.H:189
IC::Ellipsoid
Definition: Ellipsoid.H:28
Integrator::Fracture::beta
Set::Scalar beta
Definition: Fracture.H:579
Integrator::Fracture::energy_pristine_old_mf
Set::Field< Set::Scalar > energy_pristine_old_mf
energy of the pristine material for previous time step.
Definition: Fracture.H:546
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
IC::Notch
Definition: Notch.H:17
pp_queryclass
#define pp_queryclass(...)
Definition: ParmParse.H:105
Solid.H
Set::Matrix
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition: Base.H:23
Integrator::Fracture::driving_force_reference
Set::Scalar driving_force_reference
Definition: Fracture.H:563
Integrator::Fracture::ic
IC::IC * ic
Definition: Fracture.H:589
Integrator::Fracture::rhs_mf
Set::Field< Set::Scalar > rhs_mf
rhs fab for elastic solution
Definition: Fracture.H:542
Operator::Elastic::define
void define(const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const LPInfo &a_info=LPInfo(), const Vector< FabFactory< FArrayBox > const * > &a_factory={})
Definition: Elastic.cpp:26
Integrator::Fracture::driving_force_tolerance_rel
Set::Scalar driving_force_tolerance_rel
Definition: Fracture.H:566
Integrator::Fracture::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, amrex::Real, int) override
Definition: Fracture.H:482
Integrator::Fracture::nlevels
int nlevels
Definition: Fracture.H:536
Model::Solid::Linear::Isotropic
Definition: Isotropic.H:64
Sin.H
Integrator::Fracture::df_mult
Set::Scalar df_mult
mulitplier for elastic driving force.
Definition: Fracture.H:550
Crack.H
Integrator::Fracture::strain_mf
Set::Field< Set::Scalar > strain_mf
total strain field (gradient of displacement)
Definition: Fracture.H:540
Integrator::Integrator::RegisterGeneralFab
void RegisterGeneralFab(Set::Field< T > &new_fab, int ncomp, int nghost, bool evolving=true)
IC::PerturbedInterface
Definition: PerturbedInterface.H:27
Solver::Nonlocal::Newton::solve
Set::Scalar solve(const Set::Field< Set::Vector > &a_u_mf, const Set::Field< Set::Vector > &a_b_mf, Set::Field< T > &a_model_mf, Real a_tol_rel, Real a_tol_abs, const char *checkpoint_file=nullptr)
Definition: Newton.H:275
Numeric::Laplacian
AMREX_FORCE_INLINE Set::Scalar Laplacian(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const int &m, const Set::Scalar dx[AMREX_SPACEDIM])
Definition: Stencil.H:530
Integrator::Fracture::energy_mf
Set::Field< Set::Scalar > energy_mf
total elastic energy
Definition: Fracture.H:544
Integrator::Fracture::mult_df_lap
Set::Scalar mult_df_lap
Multiplier for the laplacian term.
Definition: Fracture.H:578
BC::Operator::Elastic::Constant
Definition: Constant.H:31
Integrator::Fracture::c_old_mf
Set::Field< Set::Scalar > c_old_mf
crack field at previous time step
Definition: Fracture.H:560
Integrator::Fracture::int_crack
Set::Scalar int_crack
integrated crack field over the entire domain.
Definition: Fracture.H:562
Integrator::Fracture::refinement_threshold
Set::Scalar refinement_threshold
Definition: Fracture.H:575
Constant.H
IO::ParmParse::contains
bool contains(std::string name)
Definition: ParmParse.H:151
Integrator::Fracture::elastic
struct Integrator::Fracture::@7 elastic
Integrator::Fracture::scaleModulusMax
Set::Scalar scaleModulusMax
material modulus ratio inside crack (default = 0.02).
Definition: Fracture.H:574
Integrator::Fracture::material_mf
Set::Field< Set::Scalar > material_mf
Definition: Fracture.H:584
Integrator::Fracture::energy_pristine_mf
Set::Field< Set::Scalar > energy_pristine_mf
energy of the prisitne material as if no crack is present
Definition: Fracture.H:545
Model::Interface::Crack::Constant
Definition: Constant.H:18
Integrator::Fracture::body_force
Set::Vector body_force
Definition: Fracture.H:595
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:165
Integrator::Fracture::residual_mf
Set::Field< Set::Scalar > residual_mf
residual field for solver
Definition: Fracture.H:543
Integrator
Collection of numerical integrator objects.
Definition: AllenCahn.H:39
Solver::Nonlocal::Newton::compResidual
void compResidual(Set::Field< Set::Scalar > &a_res_mf, Set::Field< Set::Scalar > &a_u_mf, Set::Field< Set::Scalar > &a_b_mf, Set::Field< T > &a_model_mf)
Definition: Newton.H:442
Operator::Operator< Grid::Node >::SetOmega
void SetOmega(Set::Scalar a_omega)
Definition: Operator.H:57
Newton.H
Integrator::Fracture::TimeStepComplete
void TimeStepComplete(amrex::Real, int) override
Definition: Fracture.H:523
Operator::Elastic::SetAverageDownCoeffs
virtual void SetAverageDownCoeffs(bool a_average_down_coeffs) override
Definition: Elastic.H:110
Operator::Elastic::Strain
void Strain(int amrlev, amrex::MultiFab &epsfab, const amrex::MultiFab &ufab, bool voigt=false) const
Definition: Elastic.cpp:403
Operator::Elastic::SetBC
void SetBC(::BC::Operator::Elastic::Elastic *a_bc)
Definition: Elastic.H:61
Operator::Elastic::Energy
void Energy(int amrlev, amrex::MultiFab &energy, const amrex::MultiFab &u, bool a_homogeneous=false)
Definition: Elastic.cpp:529
Integrator::Fracture::interval
int interval
Definition: Fracture.H:552
Integrator::Integrator::RegisterNodalFab
void RegisterNodalFab(Set::Field< Set::Scalar > &new_fab, int ncomp, int nghost, std::string name, bool writeout)
Definition: Integrator.cpp:310
Integrator::Fracture::c_mf
Set::Field< Set::Scalar > c_mf
crack field at current time step
Definition: Fracture.H:559
Integrator::Fracture::Fracture
Fracture()
Definition: Fracture.H:50
Integrator::Fracture::material
struct Integrator::Fracture::@9 material
Numeric::FieldToMatrix
AMREX_FORCE_INLINE Set::Matrix FieldToMatrix(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k)
Definition: Stencil.H:1027
Integrator::Fracture::omega
Set::Scalar omega
Definition: Fracture.H:553
Util::Warning
void Warning(std::string file, std::string func, int line, Args const &... args)
Definition: Util.H:171
BC.H
IO::ParmParse
Definition: ParmParse.H:110
Ellipsoid.H
Integrator::Fracture::stress_mf
Set::Field< Set::Scalar > stress_mf
stress field
Definition: Fracture.H:541
Numeric::Stencil
Definition: Stencil.H:38
Operator::Elastic::Stress
void Stress(int amrlev, amrex::MultiFab &sigmafab, const amrex::MultiFab &ufab, bool voigt=false, bool a_homogeneous=false)
Definition: Elastic.cpp:464
Integrator::Fracture::input_material
std::string input_material
Definition: Fracture.H:586
Integrator::Fracture::do_solve_now
bool do_solve_now
Definition: Fracture.H:551
Integrator::Fracture::cracktype
Model::Interface::Crack::Constant cracktype
type of crack. See Crack/Constant or Crack/Sin
Definition: Fracture.H:569
Integrator::Fracture::int_energy
Set::Scalar int_energy
integrated energy over the entire domain.
Definition: Fracture.H:547
IO::ParmParse::queryclass
void queryclass(std::string name, T *value, std::string file="", std::string func="", int line=-1)
Definition: ParmParse.H:443
Integrator::Integrator::integrate_variables_before_advance
bool integrate_variables_before_advance
Definition: Integrator.H:313
Integrator::Fracture::number_of_materials
int number_of_materials
Definition: Fracture.H:535
IC::Laminate
Definition: Laminate.H:16
Isotropic.H
Integrator::Fracture::brittlemodeltype
amrex::Vector< brittle_fracture_model_type_test > brittlemodeltype
Definition: Fracture.H:583
INFO
#define INFO
Definition: Util.H:20
Integrator::brittle_fracture_model_type_test
Model::Solid::Linear::Isotropic brittle_fracture_model_type_test
Definition: Fracture.H:43
Integrator::Fracture::brittlebc
BC::Operator::Elastic::Constant brittlebc
elastic BC if using brittle fracture
Definition: Fracture.H:549
Integrator::Integrator::dt
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition: Integrator.H:303
IC.H
Set::Field
Definition: Set.H:42
IO::ParmParse::queryarr
int queryarr(std::string name, std::vector< T > &value, std::string, std::string, int)
Definition: ParmParse.H:333
Integrator::Integrator::RegisterIntegratedVariable
void RegisterIntegratedVariable(Set::Scalar *integrated_variable, std::string name, bool extensive=true)
Definition: Integrator.cpp:320
Integrator::Fracture::TimeStepBegin
void TimeStepBegin(Set::Scalar, int) override
Definition: Fracture.H:221
Integrator::Fracture::ic
amrex::Vector< IC::IC * > ic
crack IC. See IC/Notch and IC/Ellipsoid
Definition: Fracture.H:571
Linear.H
Integrator::Fracture::mult_df_Gc
Set::Scalar mult_df_Gc
Multiplier for Gc/zeta term.
Definition: Fracture.H:577
Integrator::Fracture::loading
struct Integrator::Fracture::@10 loading
Util::Message
void Message(std::string file, std::string func, int line, Args const &... args)
Definition: Util.H:133
Integrator::Fracture::val
Set::Scalar val
Definition: Fracture.H:596
Integrator::Fracture::ic_type
std::string ic_type
Definition: Fracture.H:590
Integrator::Fracture::driving_force_norminf
Set::Scalar driving_force_norminf
Definition: Fracture.H:565
pp_queryarr
#define pp_queryarr(...)
Definition: ParmParse.H:103