New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 6055 – NEMO

Changeset 6055


Ignore:
Timestamp:
2015-12-15T17:19:09+01:00 (8 years ago)
Author:
gm
Message:

#1613: vvl by default : update DOC on TRA advection + change in namelist comment

Location:
branches/2015/dev_r5836_NOC3_vvl_by_default
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • branches/2015/dev_r5836_NOC3_vvl_by_default/DOC/TexFiles/Biblio/Biblio.bib

    r6040 r6055  
    763763} 
    764764 
     765@PHDTHESIS{Demange_PhD2014, 
     766  author = {J. Farge}, 
     767  title = {Sch\'{e}́mas num\'{e}́riques d’advection et de propagation d’ondes de gravit\'{e}́  
     768           dans les mod\`{e}les de circulation oc\'{e}́anique.}, 
     769  school = {Doctorat es Applied Mathematiques, Grenoble University, France}, 
     770  year = {2014}, 
     771  pages = {138pp} 
     772} 
    765773 
    766774@ARTICLE{Dobricic_al_OS07, 
     
    983991  author = {M. Farge}, 
    984992  title = {Dynamique non lineaire des ondes et des tourbillons dans les equations de Saint Venant}, 
    985   school = {Doctorat es Mathematiques, Paris VI University}, 
     993  school = {Doctorat es Mathematiques, Paris VI University, France}, 
    986994  year = {1987}, 
    987995  pages = {401pp} 
     
    16631671} 
    16641672 
     1673 
     1674@ARTICLE{Lemarie_OM2015, 
     1675  author = {F. Lemari\'{e} and L. Debreu and J. Demange and  G. Madec and J.M. Molines and M. Honnorat}, 
     1676  title = {Stability Constraints for Oceanic Numerical Models:  
     1677           Implications for the Formulation of time and space Discretizations}, 
     1678  journal = OM, 
     1679  year = {2015}, 
     1680  volume = {92},  
     1681  pages = {124--148}, 
     1682  doi = {10.1016/j.ocemod.2015.06.006}, 
     1683  url = {http://dx.doi.org/10.1016/j.ocemod.2015.06.006} 
     1684} 
     1685 
    16651686@ARTICLE{Lermusiaux2001, 
    16661687  author = {P. F. J. Lermusiaux}, 
     
    16751696  author = {M. L\'{e}vy}, 
    16761697  title = {Mod\'{e}lisation des processus biog\'{e}ochimiques en M\'{e}diterran\'{e}e 
    1677    nord-occidentale. Cycle saisonnier et variabilit\'{e} m\'{e}so\'{e}chelle}, 
     1698           nord-occidentale. Cycle saisonnier et variabilit\'{e} m\'{e}so\'{e}chelle}, 
    16781699  school = {Universit\'{e} Pierre et Marie Curie, Paris, France, 207pp}, 
    16791700  year = {1996} 
     
    18141835  year = {2010}, 
    18151836  pages = {submitted}, 
     1837} 
     1838 
     1839@ARTICLE{Lele_JCP1992, 
     1840  author = {S.K. Lele}, 
     1841  title = {Compact finite difference schemes with spectral-like resolution}, 
     1842  journal = JCP, 
     1843  year = {1992}, 
     1844  volume = {103} 
     1845  pages = {16--42} 
    18161846} 
    18171847 
  • branches/2015/dev_r5836_NOC3_vvl_by_default/DOC/TexFiles/Chapters/Chap_TRA.tex

    r6040 r6055  
    135135when those parameterisations are used (see Chap.~\ref{LDF}). 
    136136 
    137 The choice of an advection scheme is made in the \textit{\ngn{namtra\_adv}} namelist, by  
    138 setting to \textit{true} one of the logicals \textit{ln\_traadv\_xxx}. The  
    139 corresponding code can be found in the \textit{traadv\_xxx.F90} module, where  
    140 \textit{xxx} is a 3 or 4 letter acronym corresponding to each scheme.  
     137Several tracer advection scheme are proposed, namely  
     138a $2^{nd}$ or $4^{th}$ order centred schemes (CEN),  
     139a $2^{nd}$ or $4^{th}$ order Flux Corrected Transport scheme (FCT), 
     140a Monotone Upstream Scheme for Conservative Laws scheme (MUSCL), 
     141a $3^{rd}$ Upstream Biased Scheme (UBS, also often called UP3), and 
     142a Quadratic Upstream Interpolation for Convective Kinematics with  
     143Estimated Streaming Terms scheme (QUICKEST). 
     144The choice is made in the \textit{\ngn{namtra\_adv}} namelist, by  
     145setting to \textit{true} one of the logicals \textit{ln\_traadv\_xxx}.  
     146The corresponding code can be found in the \textit{traadv\_xxx.F90} module,  
     147where \textit{xxx} is a 3 or 4 letter acronym corresponding to each scheme.  
    141148By default ($i.e.$ in the reference namelist, \ngn{namelist\_ref}), all the logicals  
    142149are set to \textit{false}. If the user does not select an advection scheme  
    143 in the configuration namelist (\ngn{namelist\_cfg}), the tracers will not be advected ! 
    144  
    145 Details of the advection schemes are given below. The choice of an advection scheme  
     150in the configuration namelist (\ngn{namelist\_cfg}), the tracers will \textit{not} be advected ! 
     151 
     152Details of the advection schemes are given below. The choosing an advection scheme  
    146153is a complex matter which depends on the model physics, model resolution,  
    147 type of tracer, as well as the issue of numerical cost.  
    148  
    149 Note that  
     154type of tracer, as well as the issue of numerical cost. In particular, we note that 
    150155(1) CEN and FCT schemes require an explicit diffusion operator  
    151 while the other schemes are diffusive enough so that they do not necessarily require additional diffusion ;  
     156while the other schemes are diffusive enough so that they do not necessarily need additional diffusion ;  
    152157(2) CEN and UBS are not \textit{positive} schemes 
    153158\footnote{negative values can appear in an initially strictly positive tracer field  
     
    166171%        2nd and 4th order centred schemes 
    167172% ------------------------------------------------------------------------------------------------------------- 
    168 \subsection   [$2^{nd}$ and $4^{th}$ order centred schemes (CEN) (\np{ln\_traadv\_cen})] 
    169          {$2^{nd}$ and $4^{th}$ order centred schemes (CEN) (\np{ln\_traadv\_cen}=true)} 
     173\subsection [centred schemes (CEN) (\np{ln\_traadv\_cen})] 
     174            {centred schemes (CEN) (\np{ln\_traadv\_cen}=true)} 
    170175\label{TRA_adv_cen} 
    171176 
    172177%        2nd order centred scheme   
    173178 
    174 In the $2^{nd}$ order centred formulation (CEN2), the tracer at velocity points is  
    175 evaluated as the mean of the two neighbouring $T$-point values.  
     179The centred advection scheme (CEN) is used when \np{ln\_traadv\_cen}~=~\textit{true}.  
     180Its order ($2^{nd}$ or $4^{th}$) can be chosen independently on horizontal (iso-level)  
     181and vertical direction by setting \np{nn\_cen\_h} and \np{nn\_cen\_v} to $2$ or $4$.  
     182CEN implementation can be found in the \mdl{traadv\_cen} module. 
     183 
     184In the $2^{nd}$ order centred formulation (CEN2), the tracer at velocity points  
     185is evaluated as the mean of the two neighbouring $T$-point values.  
    176186For example, in the $i$-direction : 
    177187\begin{equation} \label{Eq_tra_adv_cen2} 
     
    185195a leapfrog scheme in conjunction with an Asselin time-filter, so $T$ in  
    186196(\ref{Eq_tra_adv_cen2}) is the \textit{now} tracer value.  
    187 CEN2 is computed in the \mdl{traadv\_cen} module. 
    188197 
    189198Note that using the CEN2, the overall tracer advection is of second  
    190199order accuracy since both (\ref{Eq_tra_adv}) and (\ref{Eq_tra_adv_cen2})  
    191 have this order of accuracy. \gmcomment{Note also that ... blah, blah} 
     200have this order of accuracy. 
    192201 
    193202%        4nd order centred scheme   
    194203 
    195 In the $4^{th}$ order formulation (CEN4), tracer values are evaluated at velocity points as  
     204In the $4^{th}$ order formulation (CEN4), tracer values are evaluated at u- and v-points as  
    196205a $4^{th}$ order interpolation, and thus depend on the four neighbouring $T$-points.  
    197206For example, in the $i$-direction: 
     
    200209=\overline{   T - \frac{1}{6}\,\delta _i \left[ \delta_{i+1/2}[T] \,\right]   }^{\,i+1/2} 
    201210\end{equation} 
     211In the vertical direction (\np{nn\_cen\_v}=$4$), a $4^{th}$ COMPACT interpolation  
     212has been prefered \citep{Demange_PhD2014}. 
     213In the COMPACT scheme, both the field and its derivative are interpolated,  
     214which leads, after a matrix inversion, spectral characteristics  
     215similar to schemes of higher order \citep{Lele_JCP1992}. 
     216  
    202217 
    203218Strictly speaking, the CEN4 scheme is not a $4^{th}$ order advection scheme  
     
    212227 
    213228A direct consequence of the pseudo-fourth order nature of the scheme is that  
    214 it is not non-diffusive, $i.e.$ the global variance of a tracer is not preserved using  
    215 CEN4. Furthermore, it must be used in conjunction with an explicit  
    216 diffusion operator to produce a sensible solution. As in CEN2 case, the time-stepping is  
    217 performed using a leapfrog scheme in conjunction with an Asselin time-filter,  
    218 so $T$ in (\ref{Eq_tra_adv_cen4}) is the \textit{now} tracer. 
    219  
    220 At a $T$-grid cell adjacent to a boundary (coastline, bottom and surface), an  
    221 additional hypothesis must be made to evaluate $\tau _u^{cen4}$. This  
    222 hypothesis usually reduces the order of the scheme. Here we choose to set  
    223 the gradient of $T$ across the boundary to zero. Alternative conditions can be  
    224 specified, such as a reduction to a second order scheme for these near boundary  
    225 grid points. 
     229it is not non-diffusive, $i.e.$ the global variance of a tracer is not preserved using CEN4.  
     230Furthermore, it must be used in conjunction with an explicit diffusion operator  
     231to produce a sensible solution.  
     232As in CEN2 case, the time-stepping is performed using a leapfrog scheme in conjunction  
     233with an Asselin time-filter, so $T$ in (\ref{Eq_tra_adv_cen4}) is the \textit{now} tracer. 
     234 
     235At a $T$-grid cell adjacent to a boundary (coastline, bottom and surface),  
     236an additional hypothesis must be made to evaluate $\tau _u^{cen4}$.  
     237This hypothesis usually reduces the order of the scheme.  
     238Here we choose to set the gradient of $T$ across the boundary to zero.  
     239Alternative conditions can be specified, such as a reduction to a second order scheme  
     240for these near boundary grid points. 
    226241 
    227242% ------------------------------------------------------------------------------------------------------------- 
    228243%        FCT scheme   
    229244% ------------------------------------------------------------------------------------------------------------- 
    230 \subsection   [$2^{nd}$ and $4^{th}$ Flux Corrected Transport schemes (FCT) (\np{ln\_traadv\_fct})] 
    231          {$2^{nd}$ and $4^{th}$ Flux Corrected Transport schemes (FCT) (\np{ln\_traadv\_fct}=true)} 
     245\subsection   [Flux Corrected Transport schemes (FCT) (\np{ln\_traadv\_fct})] 
     246         {Flux Corrected Transport schemes (FCT) (\np{ln\_traadv\_fct}=true)} 
    232247\label{TRA_adv_tvd} 
    233248 
    234 In the Flux Corrected Transport formulation, the tracer at velocity  
    235 points is evaluated using a combination of an upstream and a centred scheme.  
    236 For example, in the $i$-direction : 
     249The Flux Corrected Transport schemes (FCT) is used when \np{ln\_traadv\_fct}~=~\textit{true}.  
     250Its order ($2^{nd}$ or $4^{th}$) can be chosen independently on horizontal (iso-level)  
     251and vertical direction by setting \np{nn\_fct\_h} and \np{nn\_fct\_v} to $2$ or $4$. 
     252FCT implementation can be found in the \mdl{traadv\_fct} module. 
     253 
     254In FCT formulation, the tracer at velocity points is evaluated using a combination of  
     255an upstream and a centred scheme. For example, in the $i$-direction : 
    237256\begin{equation} \label{Eq_tra_adv_fct} 
    238257\begin{split} 
     
    246265\end{equation} 
    247266where $c_u$ is a flux limiter function taking values between 0 and 1.  
     267The FCT order is the one of the centred scheme used ($i.e.$ it depends on the setting of 
     268\np{nn\_fct\_h} and \np{nn\_fct\_v}. 
    248269There exist many ways to define $c_u$, each corresponding to a different  
    249 total variance decreasing scheme. The one chosen in \NEMO is described in  
    250 \citet{Zalesak_JCP79}. $c_u$ only departs from $1$ when the advective term  
    251 produces a local extremum in the tracer field. The resulting scheme is quite  
    252 expensive but \emph{positive}. It can be used on both active and passive tracers.  
    253 This scheme is tested and compared with MUSCL and a MPDATA scheme in \citet{Levy_al_GRL01}.  
    254 The FCT scheme is implemented in the \mdl{traadv\_fct} module. 
    255  
    256 For stability reasons (see \S\ref{STP}), 
    257 $\tau _u^{cen}$ is evaluated  in (\ref{Eq_tra_adv_fct}) using the \textit{now} tracer  
    258 while $\tau _u^{ups}$ is evaluated using the \textit{before} tracer. In other words,  
     270FCT scheme. The one chosen in \NEMO is described in \citet{Zalesak_JCP79}.  
     271$c_u$ only departs from $1$ when the advective term produces a local extremum in the tracer field.  
     272The resulting scheme is quite expensive but \emph{positive}.  
     273It can be used on both active and passive tracers.  
     274A comparison of FCT-2 with MUSCL and a MPDATA scheme can be found in \citet{Levy_al_GRL01}.  
     275 
     276An additional option has been added controlled by \np{nn\_fct\_zts}. By setting this integer to  
     277a value larger than zero, a $2^{nd}$ order FCT scheme is used on both horizontal and vertical direction,  
     278but on the latter, a split-explicit time stepping is used, with a number of sub-timestep equals 
     279to \np{nn\_fct\_zts}. This option can be useful when the size of the timestep is limited  
     280by vertical advection \citep{Lemarie_OM2015)}. Note that in this case, a similar split-explicit  
     281time stepping should be used on vertical advection of momentum to insure a better stability 
     282(see \S\ref{DYN_zad}). 
     283 
     284For stability reasons (see \S\ref{STP}), $\tau _u^{cen}$ is evaluated in (\ref{Eq_tra_adv_fct})  
     285using the \textit{now} tracer while $\tau _u^{ups}$ is evaluated using the \textit{before} tracer. In other words,  
    259286the advective part of the scheme is time stepped with a leap-frog scheme  
    260287while a forward scheme is used for the diffusive part.  
     
    267294\label{TRA_adv_mus} 
    268295 
    269 The Monotone Upstream Scheme for Conservative Laws (MUSCL) has been  
    270 implemented by \citet{Levy_al_GRL01}. In its formulation, the tracer at velocity points  
     296The Monotone Upstream Scheme for Conservative Laws (MUSCL) is used when \np{ln\_traadv\_mus}~=~\textit{true}.  
     297MUSCL implementation can be found in the \mdl{traadv\_mus} module. 
     298 
     299MUSCL has been first implemented in \NEMO by \citet{Levy_al_GRL01}. In its formulation, the tracer at velocity points  
    271300is evaluated assuming a linear tracer variation between two $T$-points  
    272301(Fig.\ref{Fig_adv_scheme}). For example, in the $i$-direction : 
     
    288317directed toward land, an upstream flux is used. This choice ensure  
    289318the \textit{positive} character of the scheme.  
     319In addition, fluxes round a grid-point where a runoff is applied can optionally be  
     320computed using upstream fluxes (\np{ln\_mus\_ups}~=~\textit{true}). 
    290321 
    291322% ------------------------------------------------------------------------------------------------------------- 
     
    296327\label{TRA_adv_ubs} 
    297328 
    298 The UBS advection scheme (also often called UP3) is an upstream-biased third order  
    299 scheme based on an upstream-biased parabolic interpolation. It is also known as  
    300 the Cell Averaged QUICK scheme (Quadratic Upstream Interpolation for Convective  
    301 Kinematics). For example, in the $i$-direction : 
     329The Upstream-Biased Scheme (UBS) is used when \np{ln\_traadv\_ubs}~=~\textit{true}.  
     330UBS implementation can be found in the \mdl{traadv\_mus} module. 
     331 
     332The UBS scheme, often called UP3, is also known as the Cell Averaged QUICK scheme  
     333(Quadratic Upstream Interpolation for Convective Kinematics). It is an upstream-biased  
     334third order scheme based on an upstream-biased parabolic interpolation.   
     335For example, in the $i$-direction : 
    302336\begin{equation} \label{Eq_tra_adv_ubs} 
    303337   \tau _u^{ubs} =\overline T ^{i+1/2}-\;\frac{1}{6} \left\{       
     
    313347 the advection scheme is similar to that reported in \cite{Farrow1995}.  
    314348It is a relatively good compromise between accuracy and smoothness.  
    315 It is not a \emph{positive} scheme, meaning that false extrema are permitted,  
     349Nevertheless the scheme is not \emph{positive}, meaning that false extrema are permitted,  
    316350but the amplitude of such are significantly reduced over the centred second  
    317 or fourth order method. Nevertheless it is not recommended that it should be  
     351or fourth order method. therefore it is not recommended that it should be  
    318352applied to a passive tracer that requires positivity.  
    319353 
    320354The intrinsic diffusion of UBS makes its use risky in the vertical direction  
    321 where the control of artificial diapycnal fluxes is of paramount importance.  
    322 Therefore the vertical flux is evaluated using either a 2nd order FCT scheme  
    323 or a 4th order COMPACT scheme (\np{nn\_cen\_v}=2 or 4). 
     355where the control of artificial diapycnal fluxes is of paramount importance \citep{Shchepetkin_McWilliams_OM05, Demange_PhD2014}.  
     356Therefore the vertical flux is evaluated using either a $2^nd$ order FCT scheme  
     357or a $4^th$ order COMPACT scheme (\np{nn\_cen\_v}=2 or 4). 
    324358 
    325359For stability reasons  (see \S\ref{STP}), 
     
    336370substitution in the \mdl{traadv\_ubs} module and obtain a QUICK scheme. 
    337371 
    338 ??? 
    339  
    340 Four different options are possible for the vertical  
    341 component used in the UBS scheme. $\tau _w^{ubs}$ can be evaluated  
    342 using either \textit{(a)} a centred $2^{nd}$ order scheme, or  \textit{(b)}  
    343 a FCT scheme, or  \textit{(c)} an interpolation based on conservative  
    344 parabolic splines following the \citet{Shchepetkin_McWilliams_OM05}  
    345 implementation of UBS in ROMS, or  \textit{(d)} a UBS. The $3^{rd}$ case  
    346 has dispersion properties similar to an eighth-order accurate conventional scheme. 
    347 The current reference version uses method (b). 
    348  
    349 ??? 
    350  
    351 Note that : 
    352  
    353 (1) When a high vertical resolution $O(1m)$ is used, the model stability can  
    354 be controlled by vertical advection (not vertical diffusion which is usually  
    355 solved using an implicit scheme). Computer time can be saved by using a  
    356 time-splitting technique on vertical advection. Such a technique has been  
    357 implemented and validated in ORCA05 with 301 levels. It is not available  
    358 in the current reference version.  
    359  
    360 (2) It is straightforward to rewrite \eqref{Eq_tra_adv_ubs} as follows: 
     372Note that it is straightforward to rewrite \eqref{Eq_tra_adv_ubs} as follows: 
    361373\begin{equation} \label{Eq_traadv_ubs2} 
    362374\tau _u^{ubs} = \tau _u^{cen4} + \frac{1}{12} \left\{   
     
    380392Thirdly, the diffusion term is in fact a biharmonic operator with an eddy  
    381393coefficient which is simply proportional to the velocity: 
    382  $A_u^{lm}= \frac{1}{12}\,{e_{1u}}^3\,|u|$. Note the current version of NEMO still uses  
    383  \eqref{Eq_tra_adv_ubs}, not \eqref{Eq_traadv_ubs2}. 
    384  %%% 
    385  \gmcomment{the change in UBS scheme has to be done} 
    386  %%% 
     394 $A_u^{lm}= \frac{1}{12}\,{e_{1u}}^3\,|u|$. Note the current version of NEMO uses  
     395the computationally more efficient formulation \eqref{Eq_tra_adv_ubs}. 
    387396 
    388397% ------------------------------------------------------------------------------------------------------------- 
     
    395404The Quadratic Upstream Interpolation for Convective Kinematics with  
    396405Estimated Streaming Terms (QUICKEST) scheme proposed by \citet{Leonard1979}  
    397 is the third order Godunov scheme. It is associated with the ULTIMATE QUICKEST  
     406is used when \np{ln\_traadv\_qck}~=~\textit{true}.  
     407QUICKEST implementation can be found in the \mdl{traadv\_mus} module. 
     408 
     409QUICKEST is the third order Godunov scheme which is associated with the ULTIMATE QUICKEST  
    398410limiter \citep{Leonard1991}. It has been implemented in NEMO by G. Reffray  
    399411(MERCATOR-ocean) and can be found in the \mdl{traadv\_qck} module. 
     
    405417This no longer guarantees the positivity of the scheme. The use of TVD in the vertical  
    406418direction (as for the UBS case) should be implemented to restore this property. 
     419 
     420%%%gmcomment   :  Cross term are missing in the current implementation.... 
    407421 
    408422 
     
    12411255\hline 
    12421256coeff.   & computer name   & S-EOS     &  description                      \\ \hline 
    1243 $a_0$       & \np{nn\_a0}     & 1.6550 $10^{-1}$ &  linear thermal expansion coeff.    \\ \hline 
    1244 $b_0$       & \np{nn\_b0}     & 7.6554 $10^{-1}$ &  linear haline  expansion coeff.    \\ \hline 
    1245 $\lambda_1$ & \np{nn\_lambda1}& 5.9520 $10^{-2}$ &  cabbeling coeff. in $T^2$          \\ \hline 
    1246 $\lambda_2$ & \np{nn\_lambda2}& 5.4914 $10^{-4}$ &  cabbeling coeff. in $S^2$       \\ \hline 
    1247 $\nu$       & \np{nn\_nu}     & 2.4341 $10^{-3}$ &  cabbeling coeff. in $T \, S$       \\ \hline 
    1248 $\mu_1$     & \np{nn\_mu1}    & 1.4970 $10^{-4}$ &  thermobaric coeff. in T         \\ \hline 
    1249 $\mu_2$     & \np{nn\_mu2}    & 1.1090 $10^{-5}$ &  thermobaric coeff. in S            \\ \hline 
     1257$a_0$       & \np{rn\_a0}     & 1.6550 $10^{-1}$ &  linear thermal expansion coeff.    \\ \hline 
     1258$b_0$       & \np{rn\_b0}     & 7.6554 $10^{-1}$ &  linear haline  expansion coeff.    \\ \hline 
     1259$\lambda_1$ & \np{rn\_lambda1}& 5.9520 $10^{-2}$ &  cabbeling coeff. in $T^2$          \\ \hline 
     1260$\lambda_2$ & \np{rn\_lambda2}& 5.4914 $10^{-4}$ &  cabbeling coeff. in $S^2$       \\ \hline 
     1261$\nu$       & \np{rn\_nu}     & 2.4341 $10^{-3}$ &  cabbeling coeff. in $T \, S$       \\ \hline 
     1262$\mu_1$     & \np{rn\_mu1}    & 1.4970 $10^{-4}$ &  thermobaric coeff. in T         \\ \hline 
     1263$\mu_2$     & \np{rn\_mu2}    & 1.1090 $10^{-5}$ &  thermobaric coeff. in S            \\ \hline 
    12501264\end{tabular} 
    12511265\caption{ \label{Tab_SEOS} 
  • branches/2015/dev_r5836_NOC3_vvl_by_default/DOC/TexFiles/Namelist/namtra_adv

    r6040 r6055  
    1313      ln_mus_ups =  .false.         !  use upstream scheme near river mouths 
    1414   ln_traadv_ubs =  .false.  !  UBS scheme 
    15       nn_ubs_v   =  2               !  =2  , vertical 2nd order FCT 
     15      nn_ubs_v   =  2               !  =2  , vertical 2nd order FCT / COMPACT 4th order 
    1616   ln_traadv_qck =  .false.  !  QUICKEST scheme 
    1717/ 
  • branches/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/CONFIG/SHARED/namelist_ref

    r6004 r6055  
    765765      ln_mus_ups =  .false.         !  use upstream scheme near river mouths 
    766766   ln_traadv_ubs =  .false.  !  UBS scheme 
    767       nn_ubs_v   =  2               !  =2  , vertical 2nd order FCT 
     767      nn_ubs_v   =  2               !  =2  , vertical 2nd order FCT / COMPACT 4th order 
    768768   ln_traadv_qck =  .false.  !  QUICKEST scheme 
    769769/ 
Note: See TracChangeset for help on using the changeset viewer.