2023-03-27 10:16:23 -07:00
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
< html xmlns = "http://www.w3.org/1999/xhtml" >
< head >
< meta http-equiv = "Content-Type" content = "text/xhtml;charset=UTF-8" / >
< meta http-equiv = "X-UA-Compatible" content = "IE=9" / >
< meta name = "generator" content = "Doxygen 1.8.17" / >
2023-04-23 12:52:57 -07:00
< meta name = "description" content = "PhasicFlow is an open-source parallel DEM (discrete element method) package for simulating granular flow. It is developed in C++ and can be exectued on both GPU (like CUDA) and CPU." >
2023-03-27 10:16:23 -07:00
< title > PhasicFlow: solvers/sphereGranFlow/sphereGranFlow.cpp Source File< / title >
< link href = "tabs.css" rel = "stylesheet" type = "text/css" / >
< script type = "text/javascript" src = "jquery.js" > < / script >
< script type = "text/javascript" src = "dynsections.js" > < / script >
< link href = "navtree.css" rel = "stylesheet" type = "text/css" / >
< script type = "text/javascript" src = "resize.js" > < / script >
< script type = "text/javascript" src = "navtreedata.js" > < / script >
< script type = "text/javascript" src = "navtree.js" > < / script >
< link href = "search/search.css" rel = "stylesheet" type = "text/css" / >
< script type = "text/javascript" src = "search/searchdata.js" > < / script >
< script type = "text/javascript" src = "search/search.js" > < / script >
< script type = "text/javascript" >
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3& dn=gpl-2.0.txt GPL-v2 */
$(document).ready(function() { init_search(); });
/* @license-end */
< / script >
2023-04-02 11:35:43 -07:00
< script type = "text/x-mathjax-config" >
MathJax.Hub.Config({
extensions: ["tex2jax.js"],
jax: ["input/TeX","output/HTML-CSS"],
});
< / script >
< script type = "text/javascript" async = "async" src = "http://cdn.mathjax.org/mathjax/latest/MathJax.js" > < / script >
2023-03-27 10:16:23 -07:00
< link href = "doxygen.css" rel = "stylesheet" type = "text/css" / >
< link href = "customdoxygen.css" rel = "stylesheet" type = "text/css" / >
< / head >
< body >
< div id = "top" > <!-- do not remove this div, it is closed by doxygen! -->
< div id = "titlearea" >
2023-04-14 10:28:41 -07:00
< table cellspacing = "0" >
2023-03-27 10:16:23 -07:00
< tbody >
2023-04-14 10:28:41 -07:00
< tr >
< td id = "projectlogo" > < a href = "https://github.com/PhasicFlow" > < img alt = "Logo" src = "phasicFlow_logo.png" > < / a > < / td >
2023-03-27 10:16:23 -07:00
< td > < div id = "MSearchBox" class = "MSearchBoxInactive" >
< span class = "left" >
< img id = "MSearchSelect" src = "search/mag_sel.png"
onmouseover="return searchBox.OnSearchSelectShow()"
onmouseout="return searchBox.OnSearchSelectHide()"
alt=""/>
< input type = "text" id = "MSearchField" value = "Search" accesskey = "S"
onfocus="searchBox.OnSearchFieldFocus(true)"
onblur="searchBox.OnSearchFieldFocus(false)"
onkeyup="searchBox.OnSearchFieldChange(event)"/>
< / span > < span class = "right" >
< a id = "MSearchClose" href = "javascript:searchBox.CloseResultsWindow()" > < img id = "MSearchCloseImg" border = "0" src = "search/close.png" alt = "" / > < / a >
< / span >
< / div >
< / td >
< / tr >
2023-04-14 10:28:41 -07:00
< tr >
< td id = "projectbrief" >
< a href = "https://https://cemf.ir" > www.cemf.ir< / a >
< / td >
< / tr >
2023-03-27 10:16:23 -07:00
< / tbody >
< / table >
< / div >
<!-- end header part -->
<!-- Generated by Doxygen 1.8.17 -->
< script type = "text/javascript" >
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3& dn=gpl-2.0.txt GPL-v2 */
var searchBox = new SearchBox("searchBox", "search",false,'Search');
/* @license-end */
< / script >
< / div > <!-- top -->
< div id = "side-nav" class = "ui-resizable side-nav-resizable" >
< div id = "nav-tree" >
< div id = "nav-tree-contents" >
< div id = "nav-sync" class = "sync" > < / div >
< / div >
< / div >
< div id = "splitbar" style = "-moz-user-select:none;"
class="ui-resizable-handle">
< / div >
< / div >
< script type = "text/javascript" >
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3& dn=gpl-2.0.txt GPL-v2 */
$(document).ready(function(){initNavTree('sphereGranFlow_8cpp_source.html',''); initResizable(); });
/* @license-end */
< / script >
< div id = "doc-content" >
<!-- window showing the filter options -->
< div id = "MSearchSelectWindow"
onmouseover="return searchBox.OnSearchSelectShow()"
onmouseout="return searchBox.OnSearchSelectHide()"
onkeydown="return searchBox.OnSearchSelectKey(event)">
< / div >
<!-- iframe showing the search results (closed by default) -->
< div id = "MSearchResultsWindow" >
< iframe src = "javascript:void(0)" frameborder = "0"
name="MSearchResults" id="MSearchResults">
< / iframe >
< / div >
< div class = "header" >
< div class = "headertitle" >
< div class = "title" > sphereGranFlow.cpp< / div > < / div >
< / div > <!-- header -->
< div class = "contents" >
< a href = "sphereGranFlow_8cpp.html" > Go to the documentation of this file.< / a > < div class = "fragment" > < div class = "line" > < a name = "l00001" > < / a > < span class = "lineno" > 1< / span >   < / div >
< div class = "line" > < a name = "l00002" > < / a > < span class = "lineno" > 2< / span >   < span class = "comment" > /*------------------------------- phasicFlow ---------------------------------< / span > < / div >
< div class = "line" > < a name = "l00003" > < / a > < span class = "lineno" > 3< / span >   < span class = "comment" > O C enter of< / span > < / div >
< div class = "line" > < a name = "l00004" > < / a > < span class = "lineno" > 4< / span >   < span class = "comment" > O O E ngineering and< / span > < / div >
< div class = "line" > < a name = "l00005" > < / a > < span class = "lineno" > 5< / span >   < span class = "comment" > O O M ultiscale modeling of< / span > < / div >
< div class = "line" > < a name = "l00006" > < / a > < span class = "lineno" > 6< / span >   < span class = "comment" > OOOOOOO F luid flow < / span > < / div >
< div class = "line" > < a name = "l00007" > < / a > < span class = "lineno" > 7< / span >   < span class = "comment" > ------------------------------------------------------------------------------< / span > < / div >
< div class = "line" > < a name = "l00008" > < / a > < span class = "lineno" > 8< / span >   < span class = "comment" > Copyright (C): www.cemf.ir< / span > < / div >
< div class = "line" > < a name = "l00009" > < / a > < span class = "lineno" > 9< / span >   < span class = "comment" > email: hamid.r.norouzi AT gmail.com< / span > < / div >
< div class = "line" > < a name = "l00010" > < / a > < span class = "lineno" > 10< / span >   < span class = "comment" > ------------------------------------------------------------------------------ < / span > < / div >
< div class = "line" > < a name = "l00011" > < / a > < span class = "lineno" > 11< / span >   < span class = "comment" > Licence:< / span > < / div >
< div class = "line" > < a name = "l00012" > < / a > < span class = "lineno" > 12< / span >   < span class = "comment" > This file is part of phasicFlow code. It is a free software for simulating < / span > < / div >
< div class = "line" > < a name = "l00013" > < / a > < span class = "lineno" > 13< / span >   < span class = "comment" > granular and multiphase flows. You can redistribute it and/or modify it under< / span > < / div >
< div class = "line" > < a name = "l00014" > < / a > < span class = "lineno" > 14< / span >   < span class = "comment" > the terms of GNU General Public License v3 or any other later versions. < / span > < / div >
< div class = "line" > < a name = "l00015" > < / a > < span class = "lineno" > 15< / span >   < span class = "comment" > < / span > < / div >
< div class = "line" > < a name = "l00016" > < / a > < span class = "lineno" > 16< / span >   < span class = "comment" > phasicFlow is distributed to help others in their research in the field of < / span > < / div >
< div class = "line" > < a name = "l00017" > < / a > < span class = "lineno" > 17< / span >   < span class = "comment" > granular and multiphase flows, but WITHOUT ANY WARRANTY; without even the< / span > < / div >
< div class = "line" > < a name = "l00018" > < / a > < span class = "lineno" > 18< / span >   < span class = "comment" > implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.< / span > < / div >
< div class = "line" > < a name = "l00019" > < / a > < span class = "lineno" > 19< / span >   < span class = "comment" > < / span > < / div >
< div class = "line" > < a name = "l00020" > < / a > < span class = "lineno" > 20< / span >   < span class = "comment" > -----------------------------------------------------------------------------*/< / span > < / div >
2025-01-10 13:02:07 +03:30
< div class = "line" > < a name = "l00033" > < / a > < span class = "lineno" > 33< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "vocabs_8hpp.html" > vocabs.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00034" > < / a > < span class = "lineno" > 34< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "phasicFlowKokkos_8hpp.html" > phasicFlowKokkos.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00035" > < / a > < span class = "lineno" > 35< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "systemControl_8hpp.html" > systemControl.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00036" > < / a > < span class = "lineno" > 36< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "commandLine_8hpp.html" > commandLine.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00037" > < / a > < span class = "lineno" > 37< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "property_8hpp.html" > property.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00038" > < / a > < span class = "lineno" > 38< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "geometry_8hpp.html" > geometry.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00039" > < / a > < span class = "lineno" > 39< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "sphereParticles_8hpp.html" > sphereParticles.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00040" > < / a > < span class = "lineno" > 40< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "interaction_8hpp.html" > interaction.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00041" > < / a > < span class = "lineno" > 41< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "Insertions_8hpp.html" > Insertions.hpp< / a > " < / span > < / div >
2023-03-27 10:16:23 -07:00
< div class = "line" > < a name = "l00042" > < / a > < span class = "lineno" > 42< / span >   < / div >
2025-01-10 13:02:07 +03:30
< div class = "line" > < a name = "l00043" > < / a > < span class = "lineno" > 43< / span >   < span class = "comment" > //#include " readControlDict.hpp" < / span > < / div >
< div class = "line" > < a name = "l00044" > < / a > < span class = "lineno" > 44< / span >   < / div >
< div class = "line" > < a name = "l00051" > < / a > < span class = "lineno" > < a class = "line" href = "sphereGranFlow_8cpp.html#a0ddf1224851353fc92bfbff6f499fa97" > 51< / a > < / span >   < span class = "keywordtype" > int< / span > < a class = "code" href = "sphereGranFlow_8cpp.html#a0ddf1224851353fc92bfbff6f499fa97" > main< / a > ( < span class = "keywordtype" > int< / span > argc, < span class = "keywordtype" > char< / span > * argv[])< / div >
< div class = "line" > < a name = "l00052" > < / a > < span class = "lineno" > 52< / span >   {< / div >
< div class = "line" > < a name = "l00053" > < / a > < span class = "lineno" > 53< / span >   < / div >
< div class = "line" > < a name = "l00054" > < / a > < span class = "lineno" > 54< / span >   < a class = "code" href = "classpFlow_1_1commandLine.html" > pFlow::commandLine< / a > cmds(< / div >
< div class = "line" > < a name = "l00055" > < / a > < span class = "lineno" > 55< / span >   < span class = "stringliteral" > " sphereGranFlow" < / span > ,< / div >
< div class = "line" > < a name = "l00056" > < / a > < span class = "lineno" > 56< / span >   < span class = "stringliteral" > " DEM solver for non-cohesive spherical particles with particle insertion " < / span > < / div >
< div class = "line" > < a name = "l00057" > < / a > < span class = "lineno" > 57< / span >   < span class = "stringliteral" > " mechanism and moving geometry" < / span > );< / div >
< div class = "line" > < a name = "l00058" > < / a > < span class = "lineno" > 58< / span >   < / div >
< div class = "line" > < a name = "l00059" > < / a > < span class = "lineno" > 59< / span >   < span class = "keywordtype" > bool< / span > isCoupling = < span class = "keyword" > false< / span > ;< / div >
< div class = "line" > < a name = "l00060" > < / a > < span class = "lineno" > 60< / span >   < / div >
< div class = "line" > < a name = "l00061" > < / a > < span class = "lineno" > 61< / span >   < span class = "keywordflow" > if< / span > (!cmds.< a class = "code" href = "classpFlow_1_1commandLine.html#af199716992f3f8bb51c89ddcca847062" > parse< / a > (argc, argv)) < span class = "keywordflow" > return< / span > 0;< / div >
< div class = "line" > < a name = "l00062" > < / a > < span class = "lineno" > 62< / span >   < / div >
< div class = "line" > < a name = "l00063" > < / a > < span class = "lineno" > 63< / span >   < span class = "comment" > // this should be palced in each main < / span > < / div >
< div class = "line" > < a name = "l00064" > < / a > < span class = "lineno" > 64< / span >   < a class = "code" href = "classpFlow_1_1processors.html#af906dbdefab1fa8e20574cfe3624a1b6" > pFlow::processors::initProcessors< / a > (argc, argv);< / div >
< div class = "line" > < a name = "l00065" > < / a > < span class = "lineno" > 65< / span >   < a class = "code" href = "namespacepFlow.html#a56ba5dd1e49ff3384363db392fb1c770" > pFlow::initialize_pFlowProcessors< / a > ();< / div >
< div class = "line" > < a name = "l00066" > < / a > < span class = "lineno" > 66< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "initialize__Control_8hpp.html" > initialize_Control.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00067" > < / a > < span class = "lineno" > 67< / span >   < / div >
< div class = "line" > < a name = "l00068" > < / a > < span class = "lineno" > 68< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "setProperty_8hpp.html" > setProperty.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00069" > < / a > < span class = "lineno" > 69< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "setSurfaceGeometry_8hpp.html" > setSurfaceGeometry.hpp< / a > " < / span > < / div >
2023-03-27 10:16:23 -07:00
< div class = "line" > < a name = "l00070" > < / a > < span class = "lineno" > 70< / span >   < / div >
2025-01-10 13:02:07 +03:30
< div class = "line" > < a name = "l00071" > < / a > < span class = "lineno" > 71< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "sphereGranFlow_2createDEMComponents_8hpp.html" > createDEMComponents.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00072" > < / a > < span class = "lineno" > 72< / span >   < / div >
< div class = "line" > < a name = "l00073" > < / a > < span class = "lineno" > 73< / span >   < a class = "code" href = "streams_8hpp.html#aeb765df06121339620670437d217fec8" > REPORT< / a > (0)< < < span class = "stringliteral" > " \nStart of time loop . . .\n" < / span > < < < a class = "code" href = "streams_8hpp.html#a1861619f2a6e102d0043a98577c8c9e8" > END_REPORT< / a > ;< / div >
2023-03-27 10:16:23 -07:00
< div class = "line" > < a name = "l00074" > < / a > < span class = "lineno" > 74< / span >   < / div >
2025-01-10 13:02:07 +03:30
< div class = "line" > < a name = "l00075" > < / a > < span class = "lineno" > 75< / span >   < span class = "keywordflow" > do< / span > < / div >
< div class = "line" > < a name = "l00076" > < / a > < span class = "lineno" > 76< / span >   {< / div >
< div class = "line" > < a name = "l00077" > < / a > < span class = "lineno" > 77< / span >   < / div >
< div class = "line" > < a name = "l00078" > < / a > < span class = "lineno" > 78< / span >   < span class = "keywordflow" > if< / span > (! < a class = "code" href = "sphereGranFlow_2createDEMComponents_8hpp.html#a84c40199c91da9a7888debd293f2d7b9" > sphInsertion< / a > .insertParticles( < / div >
< div class = "line" > < a name = "l00079" > < / a > < span class = "lineno" > 79< / span >   < a class = "code" href = "initialize__Control_8hpp.html#a4f5e4e852648762473ecd75a907417ca" > Control< / a > .time().currentIter(),< / div >
< div class = "line" > < a name = "l00080" > < / a > < span class = "lineno" > 80< / span >   < a class = "code" href = "initialize__Control_8hpp.html#a4f5e4e852648762473ecd75a907417ca" > Control< / a > .time().currentTime(),< / div >
< div class = "line" > < a name = "l00081" > < / a > < span class = "lineno" > 81< / span >   < a class = "code" href = "initialize__Control_8hpp.html#a4f5e4e852648762473ecd75a907417ca" > Control< / a > .time().dt() ) )< / div >
< div class = "line" > < a name = "l00082" > < / a > < span class = "lineno" > 82< / span >   {< / div >
< div class = "line" > < a name = "l00083" > < / a > < span class = "lineno" > 83< / span >   < a class = "code" href = "error_8hpp.html#adfe9ae1313e6913aca3f96d3eb67906e" > fatalError< / a > < < < / div >
< div class = "line" > < a name = "l00084" > < / a > < span class = "lineno" > 84< / span >   < span class = "stringliteral" > " particle insertion failed in sphereDFlow solver.\n" < / span > ;< / div >
< div class = "line" > < a name = "l00085" > < / a > < span class = "lineno" > 85< / span >   < span class = "keywordflow" > return< / span > 1;< / div >
< div class = "line" > < a name = "l00086" > < / a > < span class = "lineno" > 86< / span >   }< / div >
< div class = "line" > < a name = "l00087" > < / a > < span class = "lineno" > 87< / span >   < / div >
< div class = "line" > < a name = "l00088" > < / a > < span class = "lineno" > 88< / span >   < span class = "comment" > // set force to zero< / span > < / div >
< div class = "line" > < a name = "l00089" > < / a > < span class = "lineno" > 89< / span >   < a class = "code" href = "setSurfaceGeometry_8hpp.html#a195e279064ba2595c36f5f8d504822cb" > surfGeometry< / a > .beforeIteration();< / div >
< div class = "line" > < a name = "l00090" > < / a > < span class = "lineno" > 90< / span >   < / div >
< div class = "line" > < a name = "l00091" > < / a > < span class = "lineno" > 91< / span >   < span class = "comment" > // set force to zero, predict, particle deletion and etc. < / span > < / div >
< div class = "line" > < a name = "l00092" > < / a > < span class = "lineno" > 92< / span >   < a class = "code" href = "iterateSphereParticles_2createDEMComponents_8hpp.html#ae7f117d3f7cdd8cae2bf2913f91e4d74" > sphParticles< / a > .< a class = "code" href = "classpFlow_1_1sphereParticles.html#ada71b97666fe3f66b31690bf12633c32" > beforeIteration< / a > ();< / div >
< div class = "line" > < a name = "l00093" > < / a > < span class = "lineno" > 93< / span >   < / div >
< div class = "line" > < a name = "l00094" > < / a > < span class = "lineno" > 94< / span >   < a class = "code" href = "sphereGranFlow_2createDEMComponents_8hpp.html#affb29a66c2605b3f871b00987e41053c" > sphInteraction< / a > .beforeIteration();< / div >
< div class = "line" > < a name = "l00095" > < / a > < span class = "lineno" > 95< / span >   < / div >
< div class = "line" > < a name = "l00096" > < / a > < span class = "lineno" > 96< / span >   < a class = "code" href = "sphereGranFlow_2createDEMComponents_8hpp.html#affb29a66c2605b3f871b00987e41053c" > sphInteraction< / a > .iterate(); < / div >
< div class = "line" > < a name = "l00097" > < / a > < span class = "lineno" > 97< / span >   < / div >
< div class = "line" > < a name = "l00098" > < / a > < span class = "lineno" > 98< / span >   < a class = "code" href = "setSurfaceGeometry_8hpp.html#a195e279064ba2595c36f5f8d504822cb" > surfGeometry< / a > .iterate();< / div >
< div class = "line" > < a name = "l00099" > < / a > < span class = "lineno" > 99< / span >   < / div >
< div class = "line" > < a name = "l00100" > < / a > < span class = "lineno" > 100< / span >   < a class = "code" href = "iterateSphereParticles_2createDEMComponents_8hpp.html#ae7f117d3f7cdd8cae2bf2913f91e4d74" > sphParticles< / a > .< a class = "code" href = "classpFlow_1_1sphereParticles.html#afa767bddda52eb71cea18f755e17d559" > iterate< / a > ();< / div >
2023-03-27 10:16:23 -07:00
< div class = "line" > < a name = "l00101" > < / a > < span class = "lineno" > 101< / span >   < / div >
2025-01-10 13:02:07 +03:30
< div class = "line" > < a name = "l00102" > < / a > < span class = "lineno" > 102< / span >   < a class = "code" href = "sphereGranFlow_2createDEMComponents_8hpp.html#affb29a66c2605b3f871b00987e41053c" > sphInteraction< / a > .afterIteration();< / div >
< div class = "line" > < a name = "l00103" > < / a > < span class = "lineno" > 103< / span >   < / div >
< div class = "line" > < a name = "l00104" > < / a > < span class = "lineno" > 104< / span >   < a class = "code" href = "setSurfaceGeometry_8hpp.html#a195e279064ba2595c36f5f8d504822cb" > surfGeometry< / a > .afterIteration();< / div >
< div class = "line" > < a name = "l00105" > < / a > < span class = "lineno" > 105< / span >   < / div >
< div class = "line" > < a name = "l00106" > < / a > < span class = "lineno" > 106< / span >   < a class = "code" href = "iterateSphereParticles_2createDEMComponents_8hpp.html#ae7f117d3f7cdd8cae2bf2913f91e4d74" > sphParticles< / a > .< a class = "code" href = "classpFlow_1_1particles.html#a5ab4b6c611c3256e54f51bbfc484d58e" > afterIteration< / a > ();< / div >
< div class = "line" > < a name = "l00107" > < / a > < span class = "lineno" > 107< / span >   < / div >
2023-03-27 10:16:23 -07:00
< div class = "line" > < a name = "l00108" > < / a > < span class = "lineno" > 108< / span >   < / div >
2025-01-10 13:02:07 +03:30
< div class = "line" > < a name = "l00109" > < / a > < span class = "lineno" > 109< / span >   }< span class = "keywordflow" > while< / span > (< a class = "code" href = "initialize__Control_8hpp.html#a4f5e4e852648762473ecd75a907417ca" > Control< / a > ++);< / div >
2023-03-27 10:16:23 -07:00
< div class = "line" > < a name = "l00110" > < / a > < span class = "lineno" > 110< / span >   < / div >
2025-01-10 13:02:07 +03:30
< div class = "line" > < a name = "l00111" > < / a > < span class = "lineno" > 111< / span >   < a class = "code" href = "streams_8hpp.html#aeb765df06121339620670437d217fec8" > REPORT< / a > (0)< < < span class = "stringliteral" > " \nEnd of time loop.\n" < / span > < < < a class = "code" href = "streams_8hpp.html#a1861619f2a6e102d0043a98577c8c9e8" > END_REPORT< / a > ;< / div >
< div class = "line" > < a name = "l00112" > < / a > < span class = "lineno" > 112< / span >   < / div >
< div class = "line" > < a name = "l00113" > < / a > < span class = "lineno" > 113< / span >   < span class = "comment" > // this should be palced in each main < / span > < / div >
< div class = "line" > < a name = "l00114" > < / a > < span class = "lineno" > 114< / span >   < span class = "preprocessor" > #include " < a class = "code" href = "finalize_8hpp.html" > finalize.hpp< / a > " < / span > < / div >
< div class = "line" > < a name = "l00115" > < / a > < span class = "lineno" > 115< / span >   < a class = "code" href = "classpFlow_1_1processors.html#af07b65f79644dae88cfa7a956f6c0fca" > pFlow::processors::finalizeProcessors< / a > ();< / div >
< div class = "line" > < a name = "l00116" > < / a > < span class = "lineno" > 116< / span >   < / div >
2023-03-27 10:16:23 -07:00
< div class = "line" > < a name = "l00117" > < / a > < span class = "lineno" > 117< / span >   < / div >
2025-01-10 13:02:07 +03:30
< div class = "line" > < a name = "l00118" > < / a > < span class = "lineno" > 118< / span >   } < / div >
2023-03-27 10:16:23 -07:00
< div class = "line" > < a name = "l00119" > < / a > < span class = "lineno" > 119< / span >   < / div >
< / div > <!-- fragment --> < / div > <!-- contents -->
< / div > <!-- doc - content -->
< div class = "ttc" id = "ainitialize__Control_8hpp_html" > < div class = "ttname" > < a href = "initialize__Control_8hpp.html" > initialize_Control.hpp< / a > < / div > < / div >
2025-01-10 13:02:07 +03:30
< div class = "ttc" id = "aclasspFlow_1_1processors_html_af906dbdefab1fa8e20574cfe3624a1b6" > < div class = "ttname" > < a href = "classpFlow_1_1processors.html#af906dbdefab1fa8e20574cfe3624a1b6" > pFlow::processors::initProcessors< / a > < / div > < div class = "ttdeci" > static void initProcessors(int argc, char *argv[])< / div > < div class = "ttdoc" > Initialize MPI processors.< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "processors_8cpp_source.html#l00042" > processors.cpp:42< / a > < / div > < / div >
2023-03-27 10:16:23 -07:00
< div class = "ttc" id = "asetProperty_8hpp_html" > < div class = "ttname" > < a href = "setProperty_8hpp.html" > setProperty.hpp< / a > < / div > < / div >
2025-01-10 13:02:07 +03:30
< div class = "ttc" id = "astreams_8hpp_html_aeb765df06121339620670437d217fec8" > < div class = "ttname" > < a href = "streams_8hpp.html#aeb765df06121339620670437d217fec8" > REPORT< / a > < / div > < div class = "ttdeci" > #define REPORT(n)< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "streams_8hpp_source.html#l00039" > streams.hpp:39< / a > < / div > < / div >
2023-03-27 10:16:23 -07:00
< div class = "ttc" id = "asetSurfaceGeometry_8hpp_html" > < div class = "ttname" > < a href = "setSurfaceGeometry_8hpp.html" > setSurfaceGeometry.hpp< / a > < / div > < / div >
2025-01-10 13:02:07 +03:30
< div class = "ttc" id = "aclasspFlow_1_1processors_html_af07b65f79644dae88cfa7a956f6c0fca" > < div class = "ttname" > < a href = "classpFlow_1_1processors.html#af07b65f79644dae88cfa7a956f6c0fca" > pFlow::processors::finalizeProcessors< / a > < / div > < div class = "ttdeci" > static void finalizeProcessors()< / div > < div class = "ttdoc" > Finalize MPI processors.< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "processors_8cpp_source.html#l00088" > processors.cpp:88< / a > < / div > < / div >
2023-03-27 10:16:23 -07:00
< div class = "ttc" id = "asystemControl_8hpp_html" > < div class = "ttname" > < a href = "systemControl_8hpp.html" > systemControl.hpp< / a > < / div > < / div >
< div class = "ttc" id = "aInsertions_8hpp_html" > < div class = "ttname" > < a href = "Insertions_8hpp.html" > Insertions.hpp< / a > < / div > < / div >
2025-01-10 13:02:07 +03:30
< div class = "ttc" id = "asphereGranFlow_8cpp_html_a0ddf1224851353fc92bfbff6f499fa97" > < div class = "ttname" > < a href = "sphereGranFlow_8cpp.html#a0ddf1224851353fc92bfbff6f499fa97" > main< / a > < / div > < div class = "ttdeci" > int main(int argc, char *argv[])< / div > < div class = "ttdoc" > DEM solver for simulating granular flow of cohesion-less particles.< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "sphereGranFlow_8cpp_source.html#l00051" > sphereGranFlow.cpp:51< / a > < / div > < / div >
2023-03-31 10:54:16 -07:00
< div class = "ttc" id = "aclasspFlow_1_1commandLine_html" > < div class = "ttname" > < a href = "classpFlow_1_1commandLine.html" > pFlow::commandLine< / a > < / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "commandLine_8hpp_source.html#l00036" > commandLine.hpp:36< / a > < / div > < / div >
2025-01-10 13:02:07 +03:30
< div class = "ttc" id = "asphereGranFlow_2createDEMComponents_8hpp_html" > < div class = "ttname" > < a href = "sphereGranFlow_2createDEMComponents_8hpp.html" > createDEMComponents.hpp< / a > < / div > < / div >
< div class = "ttc" id = "ainteraction_8hpp_html" > < div class = "ttname" > < a href = "interaction_8hpp.html" > interaction.hpp< / a > < / div > < / div >
2023-03-27 10:16:23 -07:00
< div class = "ttc" id = "asetSurfaceGeometry_8hpp_html_a195e279064ba2595c36f5f8d504822cb" > < div class = "ttname" > < a href = "setSurfaceGeometry_8hpp.html#a195e279064ba2595c36f5f8d504822cb" > surfGeometry< / a > < / div > < div class = "ttdeci" > auto & surfGeometry< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "setSurfaceGeometry_8hpp_source.html#l00026" > setSurfaceGeometry.hpp:26< / a > < / div > < / div >
2025-01-10 13:02:07 +03:30
< div class = "ttc" id = "asphereGranFlow_2createDEMComponents_8hpp_html_a84c40199c91da9a7888debd293f2d7b9" > < div class = "ttname" > < a href = "sphereGranFlow_2createDEMComponents_8hpp.html#a84c40199c91da9a7888debd293f2d7b9" > sphInsertion< / a > < / div > < div class = "ttdeci" > auto sphInsertion< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "sphereGranFlow_2createDEMComponents_8hpp_source.html#l00039" > createDEMComponents.hpp:39< / a > < / div > < / div >
< div class = "ttc" id = "aphasicFlowKokkos_8hpp_html" > < div class = "ttname" > < a href = "phasicFlowKokkos_8hpp.html" > phasicFlowKokkos.hpp< / a > < / div > < / div >
2023-03-31 10:54:16 -07:00
< div class = "ttc" id = "aclasspFlow_1_1commandLine_html_af199716992f3f8bb51c89ddcca847062" > < div class = "ttname" > < a href = "classpFlow_1_1commandLine.html#af199716992f3f8bb51c89ddcca847062" > pFlow::commandLine::parse< / a > < / div > < div class = "ttdeci" > bool parse(int argc, char **argv)< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "commandLine_8cpp_source.html#l00050" > commandLine.cpp:50< / a > < / div > < / div >
2025-01-10 13:02:07 +03:30
< div class = "ttc" id = "aclasspFlow_1_1particles_html_a5ab4b6c611c3256e54f51bbfc484d58e" > < div class = "ttname" > < a href = "classpFlow_1_1particles.html#a5ab4b6c611c3256e54f51bbfc484d58e" > pFlow::particles::afterIteration< / a > < / div > < div class = "ttdeci" > bool afterIteration() override< / div > < div class = "ttdoc" > This is called in time loop, after iterate.< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "particles_8cpp_source.html#l00110" > particles.cpp:110< / a > < / div > < / div >
2023-03-27 10:16:23 -07:00
< div class = "ttc" id = "ageometry_8hpp_html" > < div class = "ttname" > < a href = "geometry_8hpp.html" > geometry.hpp< / a > < / div > < / div >
2025-01-10 13:02:07 +03:30
< div class = "ttc" id = "aclasspFlow_1_1sphereParticles_html_afa767bddda52eb71cea18f755e17d559" > < div class = "ttname" > < a href = "classpFlow_1_1sphereParticles.html#afa767bddda52eb71cea18f755e17d559" > pFlow::sphereParticles::iterate< / a > < / div > < div class = "ttdeci" > bool iterate() override< / div > < div class = "ttdoc" > iterate particles< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "sphereParticles_8cpp_source.html#l00537" > sphereParticles.cpp:537< / a > < / div > < / div >
< div class = "ttc" id = "anamespacepFlow_html_a56ba5dd1e49ff3384363db392fb1c770" > < div class = "ttname" > < a href = "namespacepFlow.html#a56ba5dd1e49ff3384363db392fb1c770" > pFlow::initialize_pFlowProcessors< / a > < / div > < div class = "ttdeci" > void initialize_pFlowProcessors()< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "pFlowProcessors_8cpp_source.html#l00010" > pFlowProcessors.cpp:10< / a > < / div > < / div >
< div class = "ttc" id = "astreams_8hpp_html_a1861619f2a6e102d0043a98577c8c9e8" > < div class = "ttname" > < a href = "streams_8hpp.html#a1861619f2a6e102d0043a98577c8c9e8" > END_REPORT< / a > < / div > < div class = "ttdeci" > #define END_REPORT< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "streams_8hpp_source.html#l00040" > streams.hpp:40< / a > < / div > < / div >
< div class = "ttc" id = "aerror_8hpp_html_adfe9ae1313e6913aca3f96d3eb67906e" > < div class = "ttname" > < a href = "error_8hpp.html#adfe9ae1313e6913aca3f96d3eb67906e" > fatalError< / a > < / div > < div class = "ttdeci" > #define fatalError< / div > < div class = "ttdoc" > Report a fatal error and exit the applicaiton.< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "error_8hpp_source.html#l00070" > error.hpp:70< / a > < / div > < / div >
< div class = "ttc" id = "asphereGranFlow_2createDEMComponents_8hpp_html_affb29a66c2605b3f871b00987e41053c" > < div class = "ttname" > < a href = "sphereGranFlow_2createDEMComponents_8hpp.html#affb29a66c2605b3f871b00987e41053c" > sphInteraction< / a > < / div > < div class = "ttdeci" > auto & sphInteraction< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "sphereGranFlow_2createDEMComponents_8hpp_source.html#l00050" > createDEMComponents.hpp:50< / a > < / div > < / div >
< div class = "ttc" id = "aclasspFlow_1_1sphereParticles_html_ada71b97666fe3f66b31690bf12633c32" > < div class = "ttname" > < a href = "classpFlow_1_1sphereParticles.html#ada71b97666fe3f66b31690bf12633c32" > pFlow::sphereParticles::beforeIteration< / a > < / div > < div class = "ttdeci" > bool beforeIteration() override< / div > < div class = "ttdoc" > before iteration step< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "sphereParticles_8cpp_source.html#l00515" > sphereParticles.cpp:515< / a > < / div > < / div >
2023-03-27 10:16:23 -07:00
< div class = "ttc" id = "asphereParticles_8hpp_html" > < div class = "ttname" > < a href = "sphereParticles_8hpp.html" > sphereParticles.hpp< / a > < / div > < / div >
< div class = "ttc" id = "aproperty_8hpp_html" > < div class = "ttname" > < a href = "property_8hpp.html" > property.hpp< / a > < / div > < / div >
2025-01-10 13:02:07 +03:30
< div class = "ttc" id = "avocabs_8hpp_html" > < div class = "ttname" > < a href = "vocabs_8hpp.html" > vocabs.hpp< / a > < / div > < / div >
< div class = "ttc" id = "aiterateSphereParticles_2createDEMComponents_8hpp_html_ae7f117d3f7cdd8cae2bf2913f91e4d74" > < div class = "ttname" > < a href = "iterateSphereParticles_2createDEMComponents_8hpp.html#ae7f117d3f7cdd8cae2bf2913f91e4d74" > sphParticles< / a > < / div > < div class = "ttdeci" > pFlow::sphereParticles sphParticles(Control, proprties)< / div > < / div >
2023-03-31 10:54:16 -07:00
< div class = "ttc" id = "acommandLine_8hpp_html" > < div class = "ttname" > < a href = "commandLine_8hpp.html" > commandLine.hpp< / a > < / div > < / div >
2023-03-27 10:16:23 -07:00
< div class = "ttc" id = "ainitialize__Control_8hpp_html_a4f5e4e852648762473ecd75a907417ca" > < div class = "ttname" > < a href = "initialize__Control_8hpp.html#a4f5e4e852648762473ecd75a907417ca" > Control< / a > < / div > < div class = "ttdeci" > auto & Control< / div > < div class = "ttdef" > < b > Definition:< / b > < a href = "initialize__Control_8hpp_source.html#l00047" > initialize_Control.hpp:47< / a > < / div > < / div >
< div class = "ttc" id = "afinalize_8hpp_html" > < div class = "ttname" > < a href = "finalize_8hpp.html" > finalize.hpp< / a > < / div > < / div >
< div id = "nav-path" class = "navpath" > <!-- id is needed for treeview function! -->
< ul >
< li class = "navelem" > < a class = "el" href = "dir_65b24c28d0f232e494405d4f9f0c5236.html" > solvers< / a > < / li > < li class = "navelem" > < a class = "el" href = "dir_ac9f1d02bd348986b458efcb0494a045.html" > sphereGranFlow< / a > < / li > < li class = "navelem" > < a class = "el" href = "sphereGranFlow_8cpp.html" > sphereGranFlow.cpp< / a > < / li >
< li class = "footer" > Generated by
< a href = "http://www.doxygen.org/index.html" >
< img class = "footer" src = "doxygen.png" alt = "doxygen" / > < / a > 1.8.17 < / li >
< / ul >
< / div >
< / body >
< / html >