REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_scale_rhs_vars.cpp
Go to the documentation of this file.
1#include <REMORA.H>
2
3using namespace amrex;
4
5void
7{
8 for (int lev=0; lev<=finest_level;lev++) {
9 MultiFab& mf_cons = *cons_new[lev];
10#ifdef _OPENMP
11#pragma omp parallel if (Gpu::notInLaunchRegion())
12#endif
13 for ( MFIter mfi(mf_cons, TilingIfNotGPU()); mfi.isValid(); ++mfi )
14 {
15 Array4<Real const> const& pm = vec_pm[lev]->array(mfi);
16 Array4<Real const> const& pn = vec_pn[lev]->array(mfi);
17 Array4<Real > const& ru = vec_ru[lev]->array(mfi);
18 Array4<Real > const& rv = vec_rv[lev]->array(mfi);
19 Array4<Real > const& ru2d = vec_ru2d[lev]->array(mfi);
20 Array4<Real > const& rv2d = vec_rv2d[lev]->array(mfi);
21
22 Box ubx = mfi.grownnodaltilebox(0,IntVect(NGROW,NGROW,0));
23 Box vbx = mfi.grownnodaltilebox(1,IntVect(NGROW,NGROW,0));
24 Box ubx2d = ubx; ubx.makeSlab(2,0);
25 Box vbx2d = vbx; vbx.makeSlab(2,0);
26
27 ParallelFor(ubx, 2, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
28 {
29 Real cff = (pm(i,j,0)+pm(i-1,j,0)) * (pn(i,j,0)+pn(i-1,j,0));
30 ru(i,j,k,n) = ru(i,j,k,n) / cff;
31 });
32
33 ParallelFor(vbx, 2, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
34 {
35 Real cff = (pm(i,j,0)+pm(i,j-1,0)) * (pn(i,j,0)+pn(i,j-1,0));
36 rv(i,j,k,n) = rv(i,j,k,n) / cff;
37 });
38
39 ParallelFor(ubx2d, 2, [=] AMREX_GPU_DEVICE (int i, int j, int , int n)
40 {
41 Real cff = (pm(i,j,0)+pm(i-1,j,0)) * (pn(i,j,0)+pn(i-1,j,0));
42 ru2d(i,j,0,n) = ru2d(i,j,0,n) / cff;
43 });
44
45 ParallelFor(vbx2d, 2, [=] AMREX_GPU_DEVICE (int i, int j, int , int n)
46 {
47 Real cff = (pm(i,j,0)+pm(i,j-1,0)) * (pn(i,j,0)+pn(i,j-1,0));
48 rv2d(i,j,0,n) = rv2d(i,j,0,n) / cff;
49 });
50 }
51 }
52}
53
54void
56{
57 for (int lev=0; lev<=finest_level;lev++) {
58 MultiFab& mf_cons = *cons_new[lev];
59#ifdef _OPENMP
60#pragma omp parallel if (Gpu::notInLaunchRegion())
61#endif
62 for ( MFIter mfi(mf_cons, TilingIfNotGPU()); mfi.isValid(); ++mfi )
63 {
64 Array4<Real const> const& pm = vec_pm[lev]->array(mfi);
65 Array4<Real const> const& pn = vec_pn[lev]->array(mfi);
66 Array4<Real > const& ru = vec_ru[lev]->array(mfi);
67 Array4<Real > const& rv = vec_rv[lev]->array(mfi);
68 Array4<Real > const& ru2d = vec_ru2d[lev]->array(mfi);
69 Array4<Real > const& rv2d = vec_rv2d[lev]->array(mfi);
70
71 Box ubx = mfi.grownnodaltilebox(0,IntVect(NGROW,NGROW,0));
72 Box vbx = mfi.grownnodaltilebox(1,IntVect(NGROW,NGROW,0));
73 Box ubx2d = ubx; ubx.makeSlab(2,0);
74 Box vbx2d = vbx; vbx.makeSlab(2,0);
75
76 ParallelFor(ubx, 2, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
77 {
78 Real cff = (pm(i,j,0)+pm(i-1,j,0)) * (pn(i,j,0)+pn(i-1,j,0));
79 ru(i,j,k,n) = ru(i,j,k,n) * cff;
80 });
81
82 ParallelFor(vbx, 2, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
83 {
84 Real cff = (pm(i,j,0)+pm(i,j-1,0)) * (pn(i,j,0)+pn(i,j-1,0));
85 rv(i,j,k,n) = rv(i,j,k,n) * cff;
86 });
87
88 ParallelFor(ubx2d, 2, [=] AMREX_GPU_DEVICE (int i, int j, int , int n)
89 {
90 Real cff = (pm(i,j,0)+pm(i-1,j,0)) * (pn(i,j,0)+pn(i-1,j,0));
91 ru2d(i,j,0,n) = ru2d(i,j,0,n) * cff;
92 });
93
94 ParallelFor(vbx2d, 2, [=] AMREX_GPU_DEVICE (int i, int j, int , int n)
95 {
96 Real cff = (pm(i,j,0)+pm(i,j-1,0)) * (pn(i,j,0)+pn(i,j-1,0));
97 rv2d(i,j,0,n) = rv2d(i,j,0,n) * cff;
98 });
99 }
100 }
101
102}
#define NGROW
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv2d
v velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:247
void scale_rhs_vars()
Scale RHS momentum variables by 1/cell area, needed before FillPatch to different levels.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:358
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive scalar
Definition REMORA.H:220
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:245
void scale_rhs_vars_inv()
Scale RHS momentum variables by cell area, needed after FillPatch to different levels.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
u velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:241
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:360
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv
v velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:243