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 ParallelFor(Box(ru), 2, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
23 {
24 Real cff = (pm(i,j,0)+pm(i-1,j,0)) * (pn(i,j,0)+pn(i-1,j,0));
25 ru(i,j,k,n) = ru(i,j,k,n) / cff;
26 });
27
28 ParallelFor(Box(rv), 2, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
29 {
30 Real cff = (pm(i,j,0)+pm(i,j-1,0)) * (pn(i,j,0)+pn(i,j-1,0));
31 rv(i,j,k,n) = rv(i,j,k,n) / cff;
32 });
33
34 ParallelFor(Box(ru2d), 2, [=] AMREX_GPU_DEVICE (int i, int j, int , int n)
35 {
36 Real cff = (pm(i,j,0)+pm(i-1,j,0)) * (pn(i,j,0)+pn(i-1,j,0));
37 ru2d(i,j,0,n) = ru2d(i,j,0,n) / cff;
38 });
39
40 ParallelFor(Box(rv2d), 2, [=] AMREX_GPU_DEVICE (int i, int j, int , int n)
41 {
42 Real cff = (pm(i,j,0)+pm(i,j-1,0)) * (pn(i,j,0)+pn(i,j-1,0));
43 rv2d(i,j,0,n) = rv2d(i,j,0,n) / cff;
44 });
45 }
46 }
47}
48
49void
51{
52 for (int lev=0; lev<=finest_level;lev++) {
53 MultiFab& mf_cons = *cons_new[lev];
54#ifdef _OPENMP
55#pragma omp parallel if (Gpu::notInLaunchRegion())
56#endif
57 for ( MFIter mfi(mf_cons, TilingIfNotGPU()); mfi.isValid(); ++mfi )
58 {
59 Array4<Real const> const& pm = vec_pm[lev]->array(mfi);
60 Array4<Real const> const& pn = vec_pn[lev]->array(mfi);
61 Array4<Real > const& ru = vec_ru[lev]->array(mfi);
62 Array4<Real > const& rv = vec_rv[lev]->array(mfi);
63 Array4<Real > const& ru2d = vec_ru2d[lev]->array(mfi);
64 Array4<Real > const& rv2d = vec_rv2d[lev]->array(mfi);
65
66 ParallelFor(Box(ru), 2, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
67 {
68 Real cff = (pm(i,j,0)+pm(i-1,j,0)) * (pn(i,j,0)+pn(i-1,j,0));
69 ru(i,j,k,n) = ru(i,j,k,n) * cff;
70 });
71
72 ParallelFor(Box(rv), 2, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
73 {
74 Real cff = (pm(i,j,0)+pm(i,j-1,0)) * (pn(i,j,0)+pn(i,j-1,0));
75 rv(i,j,k,n) = rv(i,j,k,n) * cff;
76 });
77
78 ParallelFor(Box(ru2d), 2, [=] AMREX_GPU_DEVICE (int i, int j, int , int n)
79 {
80 Real cff = (pm(i,j,0)+pm(i-1,j,0)) * (pn(i,j,0)+pn(i-1,j,0));
81 ru2d(i,j,0,n) = ru2d(i,j,0,n) * cff;
82 });
83
84 ParallelFor(Box(rv2d), 2, [=] AMREX_GPU_DEVICE (int i, int j, int , int n)
85 {
86 Real cff = (pm(i,j,0)+pm(i,j-1,0)) * (pn(i,j,0)+pn(i,j-1,0));
87 rv2d(i,j,0,n) = rv2d(i,j,0,n) * cff;
88 });
89 }
90 }
91
92}
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv2d
v velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:244
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:355
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive scalar
Definition REMORA.H:217
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:242
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:238
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:357
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv
v velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:240