List of commits:
Subject Hash Author Date (UTC)
Added lower bound for the efficiency term in buoyancyInduced kernel 37484d897975e1c3ecb20462a68bd0cd09925f9e Anurag M 2018-04-26 13:16:23
moved lookup of fields to constructor for efficiency, added spatially inhomogeneous buoyancyInduced kernel 7a444cf89c257318fabfee2618a45a96231eebcc Anurag M 2018-04-24 17:14:42
included constant coalescence kernel for testing... 21d0518c24475f61f326fbf4942345ac7bd8106b Anurag M 2018-04-23 17:56:36
modified case batchSettlerFine 9e6d0df52ab1b47d72c44f79cbf2b9f3fe9d051f Anurag M 2018-04-23 14:43:20
cleaned directory structure, fixed conflicting references to other libraries 97b6c823c552d5c690a78677cf1fbbfd9258388b Anurag M 2018-04-23 14:31:37
Implemented latest changes from OpenQBMM code 309b00b557e56d25057cf1a8cd5b30b80b32aada Anurag Misra 2018-04-05 07:17:32
updated batchSettlerFine case d475795ed3fe338ac7a8732e41ffc4e57d61c4ff Anurag Misra 2018-02-12 17:49:22
included latest QBMM updates, renamed diameterModel and alphaContactAngle entries to avoid conflict with OpenFOAM defaults d3d8fbc70cd9754d8552ae9adf930a5723e94d1c Anurag Misra 2018-02-12 17:24:17
updated test case 845e8a7326d6e5132cc6903390233b60f44a6616 Anurag Misra 2018-01-31 16:17:55
suppressed debug message, modified case for test run beecf9d98ab64fc5d5eb00e8f061a02370c9e814 Anurag Misra 2018-01-31 16:06:31
Added coalescence kernel to solver, switched on coalescence in case dictionary, cleaned up case batchSettlerFine 33502e69ba649befc30d1206f00f53ab3ae53b6d Anurag Misra 2018-01-31 14:44:39
Suppressed coalescence due to conflict, switched on dummy aggregation model in test case to test code stability e38a0560ffb1c6e68d7dc410f03a152689477078 Anurag Misra 2018-01-30 09:13:31
first commit to update code to OpenFOAM50, no coalescence kernels included 9d70a283abb4e26983a8d9711d76d22c251fbd4c Anurag Misra 2018-01-23 10:50:12
fixed alpha BC at outlet, included new case for horizontalSettler with only 2 outlets c3ac36e6699dc3d9c5b12ec02160e855d500cfa1 Anurag Misra 2017-12-13 11:57:01
added residualAlpha to phaseModel, updated case dictionaries 6f0d0800829af8e706bb67ef7bb9ef344fab1cd0 Anurag Misra 2017-12-12 00:53:41
undo previous addition, add smallnumber to rAUs for stability e5eb3353a23efd00d36c379bac1f882e74a10643 Anurag Misra 2017-12-11 23:57:15
modified pEqn.H to limit alpha to 1e-6 07f66051823747fc6acbf519d757b864f3306433 Anurag Misra 2017-12-11 23:25:12
fixed missing dictionary entries in case horizontalSettler2D 47531fc01ecde45587480e12c454d968d65aa5c6 Anurag Misra 2017-12-06 15:48:40
changed boundary conditions and mesh in case horizontalSettler2D 091228be5c25fb6c2d0f08cc5e2cdaafef28ad89 Anurag Misra 2017-12-06 15:31:09
added fvOptions compatibility in solver ddcc6c4dfc93450a9320d59855dd1955eabb30d0 Anurag Misra 2017-12-06 11:43:33
Commit 37484d897975e1c3ecb20462a68bd0cd09925f9e - Added lower bound for the efficiency term in buoyancyInduced kernel
Author: Anurag M
Author date (UTC): 2018-04-26 13:16
Committer name: Anurag M
Committer date (UTC): 2018-04-26 13:16
Parent(s): 7a444cf89c257318fabfee2618a45a96231eebcc
Signer:
Signing key:
Signing status: N
Tree: 7f3dd25fcd12138831eb2b47e0de9eaa87068353
File Lines added Lines deleted
multiphaseEulerPbeFoam/quadratureMethods/populationBalanceModels/populationBalanceSubModels/coalescenceKernels/buoyancyInduced/buoyancyInduced.C 19 12
multiphaseEulerPbeFoam/quadratureMethods/populationBalanceModels/populationBalanceSubModels/coalescenceKernels/buoyancyInducedSI/buoyancyInducedSI.C 17 22
File multiphaseEulerPbeFoam/quadratureMethods/populationBalanceModels/populationBalanceSubModels/coalescenceKernels/buoyancyInduced/buoyancyInduced.C changed (mode: 100644) (index f9a3f71..3044e6e)
... ... Foam::populationBalanceSubModels::coalescenceKernels::buoyancyInduced
59 59 Ucrit_(dict.lookup("Ucrit")), Ucrit_(dict.lookup("Ucrit")),
60 60 thisPhaseName_(dict.lookup("thisPhaseName")), thisPhaseName_(dict.lookup("thisPhaseName")),
61 61 contPhaseName_(dict.lookup("contPhaseName")), contPhaseName_(dict.lookup("contPhaseName")),
62 Udisp_(mesh_.lookupObject<volVectorField>("U")),
62 Udisp_
63 (
64 mesh_.lookupObject<volVectorField>
65 (
66 IOobject::groupName("U", thisPhaseName_)
67 )
68 ),
63 69 Ucont_ Ucont_
64 70 ( (
65 71 mesh_.lookupObject<volVectorField> mesh_.lookupObject<volVectorField>
 
... ... Foam::populationBalanceSubModels::coalescenceKernels::buoyancyInduced::Kk
101 107 const label celli const label celli
102 108 ) const ) const
103 109 { {
104
105 110 //Info << "Printing objectRegistry Names abscissa1.mesh(): " << abscissa1.mesh().objectRegistry::names() << endl; //Info << "Printing objectRegistry Names abscissa1.mesh(): " << abscissa1.mesh().objectRegistry::names() << endl;
106 111
107 112 //Info << "Printing objectRegistry Names mesh(): " << mesh_.objectRegistry::names() << endl; //Info << "Printing objectRegistry Names mesh(): " << mesh_.objectRegistry::names() << endl;
 
... ... Foam::populationBalanceSubModels::coalescenceKernels::buoyancyInduced::Kk
119 124 dimensionedScalar Uslip_ dimensionedScalar Uslip_
120 125 ( (
121 126 "Uslip_", "Uslip_",
122 dimensionSet(0,1,-1,0,0,0,0),
123 mag(Udisp_[celli]-Ucont_[celli])
127 dimensionSet(0, 1, -1, 0, 0),
128 mag(Udisp_[celli] - Ucont_[celli])
124 129 ); );
125 130
126 // Moved Ucrit to populationBalancedictionary
127 //dimensionedScalar Ucrit_("Ucrit_", dimensionSet(0,1,-1,0,0,0,0), 0.001);
128
129 131 // Read sauter diameter // Read sauter diameter
130 //dimensionedScalar diaMin_("diaMin_", dimensionSet(0,1,0,0,0,0,0), 1.0e-06);
131 //dimensionedScalar diaMax_("diaMax_", dimensionSet(0,1,0,0,0,0,0), 1.0e-01);
132 //dimensionedScalar smallMoment2("smallMoment2", dimensionSet(0,-1,0,0,0,0,0), 1.0e-06);
132 //dimensionedScalar diaMin_("diaMin_", dimensionSet(0, 1, 0, 0, 0), 1.0e-06);
133 //dimensionedScalar diaMax_("diaMax_", dimensionSet(0, 1, 0, 0, 0), 1.0e-01);
134 //dimensionedScalar smallMoment2("smallMoment2", dimensionSet(0, -1, 0, 0, 0), 1.0e-06);
133 135
134 136 scalar diaMin_(1.0e-06); scalar diaMin_(1.0e-06);
135 137 scalar diaMax_(1.0e-01); scalar diaMax_(1.0e-01);
 
... ... Foam::populationBalanceSubModels::coalescenceKernels::buoyancyInduced::Kk
138 140 dimensionedScalar dSauter_ dimensionedScalar dSauter_
139 141 ( (
140 142 "dSauter_", "dSauter_",
141 dimensionSet(0,1,0,0,0,0,0),
143 dimensionSet(0, 1, 0, 0, 0),
142 144 min min
143 145 ( (
144 146 max max
145 147 ( (
146 moment3_[celli]/max(moment2_[celli],smallMoment2), diaMin_
148 moment3_[celli]/max(moment2_[celli], smallMoment2), diaMin_
147 149 ), ),
148 150 diaMax_ diaMax_
149 151 ) )
 
... ... Foam::populationBalanceSubModels::coalescenceKernels::buoyancyInduced::Kk
155 157 mag(Uslip_.value()*(abscissa2-abscissa1)/dSauter_.value()) mag(Uslip_.value()*(abscissa2-abscissa1)/dSauter_.value())
156 158 ); );
157 159
160 //if (std::fabs(Urelative_ - 0.00025) > SMALL)
161 // Urelative_ = 0.00025;
162
158 163 //- Evaluate efficiency term //- Evaluate efficiency term
159 164 scalar lambda_(Foam::exp(-(Urelative_/Ucrit_.value()))); scalar lambda_(Foam::exp(-(Urelative_/Ucrit_.value())));
165 lambda_ = max(lambda_, 1e-10);
166
160 167 /* /*
161 168 if (min(lambda_) < smallAbs) if (min(lambda_) < smallAbs)
162 169 { {
File multiphaseEulerPbeFoam/quadratureMethods/populationBalanceModels/populationBalanceSubModels/coalescenceKernels/buoyancyInducedSI/buoyancyInducedSI.C changed (mode: 100644) (index ead4365..1bab1f1)
... ... Foam::populationBalanceSubModels::coalescenceKernels::buoyancyInducedSI
59 59 Ucrit_(dict.lookup("Ucrit")), Ucrit_(dict.lookup("Ucrit")),
60 60 thisPhaseName_(dict.lookup("thisPhaseName")), thisPhaseName_(dict.lookup("thisPhaseName")),
61 61 contPhaseName_(dict.lookup("contPhaseName")), contPhaseName_(dict.lookup("contPhaseName")),
62 Udisp_(mesh_.lookupObject<volVectorField>("U")),
62 //Udisp_(mesh_.lookupObject<volVectorField>("U")),
63 Udisp_
64 (
65 mesh_.lookupObject<volVectorField>
66 (
67 IOobject::groupName("U", thisPhaseName_)
68 )
69 ),
63 70 Ucont_ Ucont_
64 71 ( (
65 72 mesh_.lookupObject<volVectorField> mesh_.lookupObject<volVectorField>
 
... ... Foam::populationBalanceSubModels::coalescenceKernels::buoyancyInducedSI::Kk
105 112 const label celli const label celli
106 113 ) const ) const
107 114 { {
108
115 //- Evaluate relative velocity
109 116 dimensionedScalar Uslip_ dimensionedScalar Uslip_
110 117 ( (
111 118 "Uslip_", "Uslip_",
112 dimensionSet(0,1,-1,0,0,0,0),
113 mag(Udisp_[celli]-Ucont_[celli])
119 dimensionSet(0, 1, -1, 0, 0),
120 mag(Udisp_[celli] - Ucont_[celli])
114 121 ); );
115 // Moved Ucrit to populationBalancedictionary
116 //dimensionedScalar Ucrit_("Ucrit_", dimensionSet(0,1,-1,0,0,0,0), 0.001);
117
118 // Read sauter diameter
119 //dimensionedScalar diaMin_("diaMin_", dimensionSet(0,1,0,0,0,0,0), 1.0e-06);
120 //dimensionedScalar diaMax_("diaMax_", dimensionSet(0,1,0,0,0,0,0), 1.0e-01);
121 //dimensionedScalar smallMoment2("smallMoment2", dimensionSet(0,-1,0,0,0,0,0), 1.0e-06);
122
122 123 scalar diaMin_(1.0e-06); scalar diaMin_(1.0e-06);
123 124 scalar diaMax_(1.0e-01); scalar diaMax_(1.0e-01);
124 125 scalar smallMoment2(1.0e-06); scalar smallMoment2(1.0e-06);
 
... ... Foam::populationBalanceSubModels::coalescenceKernels::buoyancyInducedSI::Kk
126 127 dimensionedScalar dSauter_ dimensionedScalar dSauter_
127 128 ( (
128 129 "dSauter_", "dSauter_",
129 dimensionSet(0,1,0,0,0,0,0),
130 dimensionSet(0, 1, 0, 0, 0),
130 131 min min
131 132 ( (
132 133 max max
133 134 ( (
134 moment3_[celli]/max(moment2_[celli],smallMoment2), diaMin_
135 moment3_[celli]/max(moment2_[celli], smallMoment2), diaMin_
135 136 ), ),
136 137 diaMax_ diaMax_
137 138 ) )
 
... ... Foam::populationBalanceSubModels::coalescenceKernels::buoyancyInducedSI::Kk
145 146
146 147 //- Evaluate efficiency term //- Evaluate efficiency term
147 148 scalar lambda_(Foam::exp(-(Urelative_/Ucrit_.value()))); scalar lambda_(Foam::exp(-(Urelative_/Ucrit_.value())));
148 /*
149 if (min(lambda_) < smallAbs)
150 {
151 Info << "Found negative values: " << min(lambda_) << endl;
152 }
153 */
154
149 lambda_ = max(lambda_, 1e-10);
150
155 151 scalar coalescenceK scalar coalescenceK
156 152 = Ck_.value() = Ck_.value()
157 153 *(0.25*Foam::constant::mathematical::pi) *(0.25*Foam::constant::mathematical::pi)
158 154 *sqr(abscissa1 + abscissa2)*Urelative_ *sqr(abscissa1 + abscissa2)*Urelative_
159 155 *(ifPorous_[celli]*1.0 + (1.0 - ifPorous_[celli])*lambda_); *(ifPorous_[celli]*1.0 + (1.0 - ifPorous_[celli])*lambda_);
160
161
156
162 157 return coalescenceK; return coalescenceK;
163 158 } }
164 159
Hints:
Before first commit, do not forget to setup your git environment:
git config --global user.name "your_name_here"
git config --global user.email "your@email_here"

Clone this repository using HTTP(S):
git clone https://rocketgit.com/user/velagala/MultiPhaseQBMM

Clone this repository using ssh (do not forget to upload a key first):
git clone ssh://rocketgit@ssh.rocketgit.com/user/velagala/MultiPhaseQBMM

Clone this repository using git:
git clone git://git.rocketgit.com/user/velagala/MultiPhaseQBMM

You are allowed to anonymously push to this repository.
This means that your pushed commits will automatically be transformed into a merge request:
... clone the repository ...
... make some changes and some commits ...
git push origin main