1 #ifndef SOLVER_NONLOCAL_NEWTON
2 #define SOLVER_NONLOCAL_NEWTON
8 #include "Numeric/Stencil.H"
61 amrex::Box domain(
linop->
Geom(lev).Domain());
62 domain.convert(amrex::IntVect::TheNodeVector());
65 for (MFIter mfi(*a_model_mf[lev],
false); mfi.isValid(); ++mfi)
67 amrex::Box bx = mfi.grownnodaltilebox();
70 amrex::Array4<const T>
const& model = a_model_mf[lev]->array(mfi);
71 amrex::Array4<const Set::Scalar>
const& u = a_u_mf[lev]->array(mfi);
72 amrex::Array4<Set::Matrix>
const& dw = a_dw_mf[lev]->array(mfi);
73 amrex::Array4<Set::Matrix4<AMREX_SPACEDIM, T::sym>>
const& ddw = a_ddw_mf[lev]->array(mfi);
76 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
83 dw(i, j, k) = model(i, j, k).DW(
gradu);
84 ddw(i, j, k) = model(i, j, k).DDW(
gradu);
89 dw(i, j, k) = model(i, j, k).DW(eps);
90 ddw(i, j, k) = model(i, j, k).DDW(eps);
95 dw(i, j, k) = model(i, j, k).DW(
F);
96 ddw(i, j, k) = model(i, j, k).DDW(
F);
109 amrex::Box domain(
linop->
Geom(lev).Domain());
110 domain.convert(amrex::IntVect::TheNodeVector());
112 const amrex::Dim3 lo = amrex::lbound(domain), hi = amrex::ubound(domain);
113 for (MFIter mfi(*a_model_mf[lev],
false); mfi.isValid(); ++mfi)
115 amrex::Box bx = mfi.nodaltilebox();
117 amrex::Array4<const Set::Scalar>
const& u = a_u_mf[lev]->array(mfi);
118 amrex::Array4<const Set::Scalar>
const& b = a_b_mf[lev]->array(mfi);
119 amrex::Array4<const Set::Matrix>
const& dw = a_dw_mf[lev]->array(mfi);
120 amrex::Array4<Set::Scalar>
const& rhs = a_rhs_mf[lev]->array(mfi);
121 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
125 if (AMREX_D_TERM(i == lo.x || i == hi.x, || j == lo.y || j == hi.y, || k == lo.z || k == hi.z))
128 Set::Vector U(
AMREX_D_DECL(u(i, j, k, 0), u(i, j, k, 1), u(i, j, k, 2)));
130 for (
int d = 0; d < AMREX_SPACEDIM; d++) rhs(i, j, k, d) = b(i, j, k, d) - ret(d);
135 for (
int p = 0; p < AMREX_SPACEDIM; p++) rhs(i, j, k, p) = b(i, j, k, p) - divdw(p);
153 amrex::Box domain(
linop->
Geom(lev).Domain());
154 domain.convert(amrex::IntVect::TheNodeVector());
157 for (MFIter mfi(*a_model_mf[lev],
false); mfi.isValid(); ++mfi)
159 amrex::Box bx = mfi.grownnodaltilebox();
162 amrex::Array4<const T>
const& model = a_model_mf[lev]->array(mfi);
163 amrex::Array4<const Set::Vector>
const& u = a_u_mf[lev]->array(mfi);
164 amrex::Array4<Set::Matrix>
const& dw = a_dw_mf[lev]->array(mfi);
165 amrex::Array4<Set::Matrix4<AMREX_SPACEDIM, T::sym>>
const& ddw = a_ddw_mf[lev]->array(mfi);
168 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
179 kinvar =
gradu + Set::Matrix::Identity();
181 dw(i, j, k) = model(i, j, k).DW(kinvar);
182 ddw(i, j, k) = model(i, j, k).DDW(kinvar);
199 amrex::Box domain(
linop->
Geom(lev).Domain());
200 domain.convert(amrex::IntVect::TheNodeVector());
202 const amrex::Dim3 lo = amrex::lbound(domain), hi = amrex::ubound(domain);
203 for (MFIter mfi(*a_model_mf[lev],
false); mfi.isValid(); ++mfi)
205 amrex::Box bx = mfi.grownnodaltilebox() & domain;
206 amrex::Array4<const Set::Vector>
const& u = a_u_mf[lev]->array(mfi);
207 amrex::Array4<const Set::Vector>
const& b = a_b_mf[lev]->array(mfi);
208 amrex::Array4<const Set::Matrix>
const& dw = a_dw_mf[lev]->array(mfi);
209 amrex::Array4<Set::Scalar>
const& rhs = a_rhs_mf[lev]->array(mfi);
212 if (
m_psi ==
nullptr)
214 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
218 if (AMREX_D_TERM(i == lo.x || i == hi.x, || j == lo.y || j == hi.y, || k == lo.z || k == hi.z))
221 Set::Vector ret =
m_elastic->
GetBC()(u(i, j, k),
gradu, dw(i, j, k), i, j, k, bx);
222 for (
int d = 0; d < AMREX_SPACEDIM; d++) rhs(i, j, k, d) = b(i, j, k)(d) - ret(d);
227 for (
int d = 0; d < AMREX_SPACEDIM; d++) rhs(i, j, k, d) = b(i, j, k)(d) - divdw(d);
236 const amrex::Dim3 boxlo = amrex::lbound(bx), boxhi = amrex::ubound(bx);
237 amrex::Array4<const Set::Scalar>
const& psi = (*m_psi)[lev]->array(mfi);
238 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
242 if (AMREX_D_TERM(i == lo.x || i == hi.x, || j == lo.y || j == hi.y, || k == lo.z || k == hi.z))
246 if (AMREX_D_TERM(i > boxlo.x && i<boxhi.x, && j>boxlo.y && j<boxhi.y, && k>boxlo.z && k < boxhi.z))
250 Set::Vector ret =
m_elastic->
GetBC()(u(i, j, k),
gradu, dw(i, j, k) * psiavg, i, j, k, bx);
251 for (
int d = 0; d < AMREX_SPACEDIM; d++) rhs(i, j, k, d) = b(i, j, k)(d) - ret(d);
256 if (AMREX_D_TERM(i > boxlo.x && i<boxhi.x, && j>boxlo.y && j<boxhi.y, && k>boxlo.z && k < boxhi.z))
260 divdw += dw(i, j, k) * gradpsi;
262 for (
int d = 0; d < AMREX_SPACEDIM; d++) rhs(i, j, k, d) = b(i, j, k)(d) - divdw(d);
278 Real a_tol_rel, Real a_tol_abs,
const char* checkpoint_file =
nullptr)
290 dsol_mf.
Define(lev, a_u_mf[lev]->boxArray(),
291 a_u_mf[lev]->DistributionMap(),
293 a_u_mf[lev]->nGrow());
294 dw_mf.
Define(lev, a_b_mf[lev]->boxArray(),
295 a_b_mf[lev]->DistributionMap(),
297 a_b_mf[lev]->nGrow());
298 ddw_mf.
Define(lev, a_b_mf[lev]->boxArray(),
299 a_b_mf[lev]->DistributionMap(),
301 a_b_mf[lev]->nGrow());
302 rhs_mf.
Define(lev, a_b_mf[lev]->boxArray(),
303 a_b_mf[lev]->DistributionMap(),
305 a_b_mf[lev]->nGrow());
307 dsol_mf[lev]->setVal(0.0);
308 dw_mf[lev]->setVal(Set::Matrix::Zero());
311 a_b_mf.
Copy(lev, *rhs_mf[lev], 0, 2);
315 for (
int nriter = 0; nriter <
m_nriters; nriter++)
327 for (
int lev = 0; lev < dsol_mf.size(); ++lev)
329 for (
int comp = 0; comp < AMREX_SPACEDIM; comp++)
331 Set::Scalar tmpcornorm = dsol_mf[lev]->norm0(comp, 0);
332 if (tmpcornorm > cornorm) cornorm = tmpcornorm;
340 if (solnorm == 0) relnorm = cornorm;
341 else relnorm = cornorm / solnorm;
344 for (
int lev = 0; lev < dsol_mf.size(); ++lev)
345 a_u_mf.
AddFrom(lev, *dsol_mf[lev], 0, 2);
359 Real a_tol_rel, Real a_tol_abs,
const char* checkpoint_file =
nullptr)
371 dsol_mf.
Define(lev, a_u_mf[lev]->boxArray(),
372 a_u_mf[lev]->DistributionMap(),
373 a_u_mf[lev]->nComp(),
374 a_u_mf[lev]->nGrow());
375 dw_mf.
Define(lev, a_b_mf[lev]->boxArray(),
376 a_b_mf[lev]->DistributionMap(),
378 a_b_mf[lev]->nGrow());
379 ddw_mf.
Define(lev, a_b_mf[lev]->boxArray(),
380 a_b_mf[lev]->DistributionMap(),
382 a_b_mf[lev]->nGrow());
383 rhs_mf.
Define(lev, a_b_mf[lev]->boxArray(),
384 a_b_mf[lev]->DistributionMap(),
385 a_b_mf[lev]->nComp(),
386 a_b_mf[lev]->nGrow());
388 dsol_mf[lev]->setVal(0.0);
389 dw_mf[lev]->setVal(Set::Matrix::Zero());
392 amrex::MultiFab::Copy(*rhs_mf[lev], *a_b_mf[lev], 0, 0, AMREX_SPACEDIM, 2);
395 for (
int nriter = 0; nriter <
m_nriters; nriter++)
408 for (
int lev = 0; lev < dsol_mf.size(); ++lev)
410 for (
int comp = 0; comp < AMREX_SPACEDIM; comp++)
412 Set::Scalar tmpcornorm = dsol_mf[lev]->norm0(comp, 0);
413 if (tmpcornorm > cornorm) cornorm = tmpcornorm;
415 Set::Scalar tmpsolnorm = a_u_mf[lev]->norm0(comp, 0);
416 if (tmpsolnorm > solnorm) solnorm = tmpsolnorm;
421 if (solnorm == 0) relnorm = cornorm;
422 else relnorm = cornorm / solnorm;
425 for (
int lev = 0; lev < dsol_mf.size(); ++lev)
426 amrex::MultiFab::Add(*a_u_mf[lev], *dsol_mf[lev], 0, 0, AMREX_SPACEDIM, 2);
449 dw_mf.resize(a_u_mf.size());
450 ddw_mf.resize(a_u_mf.size());
451 for (
int lev = 0; lev < a_u_mf.size(); lev++)
453 dw_mf.
Define(lev, a_b_mf[lev]->boxArray(),
454 a_b_mf[lev]->DistributionMap(),
455 1, a_b_mf[lev]->nGrow());
456 ddw_mf.
Define(lev, a_b_mf[lev]->boxArray(),
457 a_b_mf[lev]->DistributionMap(),
458 1, a_b_mf[lev]->nGrow());
459 dw_mf[lev]->setVal(Set::Matrix::Zero());
473 dw_mf.resize(a_u_mf.size());
474 ddw_mf.resize(a_u_mf.size());
475 res_mf.resize(a_u_mf.size());
476 for (
int lev = 0; lev < a_u_mf.size(); lev++)
478 dw_mf.
Define(lev, a_b_mf[lev]->boxArray(),
479 a_b_mf[lev]->DistributionMap(),
480 1, a_b_mf[lev]->nGrow());
481 ddw_mf.
Define(lev, a_b_mf[lev]->boxArray(),
482 a_b_mf[lev]->DistributionMap(),
483 1, a_b_mf[lev]->nGrow());
484 res_mf.
Define(lev, a_b_mf[lev]->boxArray(),
485 a_b_mf[lev]->DistributionMap(),
486 AMREX_SPACEDIM, a_b_mf[lev]->nGrow());
487 dw_mf[lev]->setVal(Set::Matrix::Zero());
492 for (
int lev = 0; lev < a_res_mf.size(); ++lev)
495 a_res_mf.
CopyFrom(lev, *res_mf[lev], 0, 2);
508 res_mf.Define(lev, a_b_mf[lev]->boxArray(), a_b_mf[lev]->DistributionMap(), AMREX_SPACEDIM, a_b_mf[lev]->nGrow());
509 sol_mf.Define(lev, a_b_mf[lev]->boxArray(), a_b_mf[lev]->DistributionMap(), AMREX_SPACEDIM, a_b_mf[lev]->nGrow());
510 rhs_mf.
Define(lev, a_b_mf[lev]->boxArray(), a_b_mf[lev]->DistributionMap(), AMREX_SPACEDIM, a_b_mf[lev]->nGrow());
512 a_u_mf.
Copy(lev, *sol_mf[lev], 0, 2);
513 a_b_mf.
Copy(lev, *rhs_mf[lev], 0, 2);
516 mlmg->compResidual(amrex::GetVecOfPtrs(res_mf), amrex::GetVecOfPtrs(sol_mf), amrex::GetVecOfConstPtrs(rhs_mf));
520 a_res_mf.
CopyFrom(lev, *res_mf[lev], 0, 2);
529 for (
int lev = 0; lev < a_u_mf.size(); lev++)
531 BL_PROFILE(
"Solver::Nonlocal::Newton::DW()");
533 const amrex::Real* DX =
linop->
Geom(lev).CellSize();
534 amrex::Box domain(
linop->
Geom(lev).Domain());
535 domain.convert(amrex::IntVect::TheNodeVector());
537 for (MFIter mfi(*a_u_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
539 const Box& bx = mfi.tilebox();
540 amrex::Array4<T>
const& C = a_model_mf[lev]->array(mfi);
541 amrex::Array4<amrex::Real>
const& w = a_w_mf[lev]->array(mfi);
542 amrex::Array4<const amrex::Real>
const& u = a_u_mf[lev]->array(mfi);
543 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
550 for (
int p = 0; p < AMREX_SPACEDIM; p++)
552 AMREX_D_TERM(
gradu(p, 0) = (
Numeric::Stencil<Set::Scalar, 1, 0, 0>::D(u, i, j, k, p, DX, sten));,
553 gradu(p, 1) = (
Numeric::Stencil<Set::Scalar, 0, 1, 0>::D(u, i, j, k, p, DX, sten));,
554 gradu(p, 2) = (
Numeric::Stencil<Set::Scalar, 0, 0, 1>::D(u, i, j, k, p, DX, sten)););
558 w(i, j, k) = C(i, j, k).W(
gradu);
560 w(i, j, k) = C(i, j, k).W(0.5 * (
gradu +
gradu.transpose()));
562 w(i, j, k) = C(i, j, k).W(
gradu + Set::Matrix::Identity());
572 for (
int lev = 0; lev < a_u_mf.size(); lev++)
574 BL_PROFILE(
"Solver::Nonlocal::Newton::DW()");
576 const amrex::Real* DX =
linop->
Geom(lev).CellSize();
577 amrex::Box domain(
linop->
Geom(lev).Domain());
578 domain.convert(amrex::IntVect::TheNodeVector());
580 for (MFIter mfi(*a_u_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
582 const Box& bx = mfi.tilebox();
583 amrex::Array4<T>
const& C = a_model_mf[lev]->array(mfi);
584 amrex::Array4<amrex::Real>
const& dw = a_dw_mf[lev]->array(mfi);
585 amrex::Array4<const amrex::Real>
const& u = a_u_mf[lev]->array(mfi);
586 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
593 for (
int p = 0; p < AMREX_SPACEDIM; p++)
595 AMREX_D_TERM(
gradu(p, 0) = (
Numeric::Stencil<Set::Scalar, 1, 0, 0>::D(u, i, j, k, p, DX, sten));,
596 gradu(p, 1) = (
Numeric::Stencil<Set::Scalar, 0, 1, 0>::D(u, i, j, k, p, DX, sten));,
597 gradu(p, 2) = (
Numeric::Stencil<Set::Scalar, 0, 0, 1>::D(u, i, j, k, p, DX, sten)););
603 sig = C(i, j, k).DW(
gradu);
605 sig = C(i, j, k).DW(0.5 * (
gradu +
gradu.transpose()));
607 sig = C(i, j, k).DW(
gradu + Set::Matrix::Identity());
611 AMREX_D_PICK(dw(i, j, k, 0) = sig(0, 0);
613 dw(i, j, k, 0) = sig(0, 0); dw(i, j, k, 1) = sig(0, 1);
614 dw(i, j, k, 2) = sig(1, 0); dw(i, j, k, 3) = sig(1, 1);
616 dw(i, j, k, 0) = sig(0, 0); dw(i, j, k, 1) = sig(0, 1); dw(i, j, k, 2) = sig(0, 2);
617 dw(i, j, k, 3) = sig(1, 0); dw(i, j, k, 4) = sig(1, 1); dw(i, j, k, 5) = sig(1, 2);
618 dw(i, j, k, 6) = sig(2, 0); dw(i, j, k, 7) = sig(2, 1); dw(i, j, k, 8) = sig(2, 2););