ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_Orbit.cpp File Reference
#include <ERF_Orbit.H>
Include dependency graph for ERF_Orbit.cpp:

Functions

void zenith (int &calday, amrex::MultiFab *clat, amrex::MultiFab *clon, const amrex::Vector< int > &rank_offsets, real1d &coszrs, int &ncol, const Real &eccen, const Real &mvelpp, const Real &lambm0, const Real &obliqr, amrex::Real uniform_angle)
 

Function Documentation

◆ zenith()

void zenith ( int &  calday,
amrex::MultiFab *  clat,
amrex::MultiFab *  clon,
const amrex::Vector< int > &  rank_offsets,
real1d &  coszrs,
int &  ncol,
const Real &  eccen,
const Real &  mvelpp,
const Real &  lambm0,
const Real &  obliqr,
amrex::Real  uniform_angle 
)
17 {
18  Real delta; // Solar declination angle in radians
19  Real eccf; // Earth orbit eccentricity factor
20 
21  // Populate delta & eccf
22  shr_orb_decl(calday, eccen, mvelpp, lambm0, obliqr, delta, eccf);
23 
24  // If we have a valid pointer, go through the whole machinery
25  if (clat) {
26  for (MFIter mfi(*clat); mfi.isValid(); ++mfi) {
27  const auto& tbx = mfi.tilebox();
28  auto nx = tbx.length(0);
29 
30  auto lat_array = clat->array(mfi);
31  auto lon_array = clon->array(mfi);
32 
33  const int offset = rank_offsets[mfi.LocalIndex()];
34 
35  // NOTE: lat/lon are 2D multifabs!
36  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int /*k*/)
37  {
38  auto icol = (j-tbx.smallEnd(1))*nx + (i-tbx.smallEnd(0)) + 1 + offset;
39  coszrs(icol) = shr_orb_cosz(calday, lat_array(i,j,0), lon_array(i,j,0), delta, uniform_angle);
40  });
41  }
42  }
43  // Use constant value near center of USA or the uniform angle
44  else {
45  Real cons_lat = 40.0;
46  Real cons_lon = 100.0;
47  Real val = shr_orb_cosz(calday, cons_lat, cons_lon, delta, uniform_angle);
48  yakl::memset(coszrs, val);
49  }
50 }
AMREX_GPU_HOST AMREX_FORCE_INLINE amrex::Real shr_orb_cosz(const amrex::Real &jday, const amrex::Real &lat, const amrex::Real &lon, const amrex::Real &declin, amrex::Real uniform_angle)
Definition: ERF_Orbit.H:26
AMREX_GPU_HOST AMREX_FORCE_INLINE void shr_orb_decl(const amrex::Real &calday, const amrex::Real &eccen, const amrex::Real &mvelpp, const amrex::Real &lambm0, const amrex::Real &obliqr, amrex::Real &delta, amrex::Real &eccf)
Definition: ERF_Orbit.H:52
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28

Referenced by Radiation::run().

Here is the call graph for this function:
Here is the caller graph for this function: