Interpolation scheme for electric current / heat flux
Functionality to add/problem to solve
An electric current is calculated as
The gradient of the electric potential can not be calculated using the Gauss theorem with LINEAR interpolation of the potential . A special interpolation scheme is necessary. As described below, the potential at the face needs to be weighted by the conductivity of both cells.
Target audience
The mentioned interpolation scheme is necessary for computing the gradient of
- electric potential (i.e. electric current)
- temperature (i.e. heat flux)
What does success look like, and how can we measure that?
I use the mentioned interpolation scheme for several years. It is already validated and published.
Links / references
Weber, N.; Beckstein, P.; Galindo, V.; Starace, M.; Weier, T.: Electro-vortex flow simulation using coupled meshes, Computers and Fluids 168(2018) 101-109
It is Eqn. 16/17 in the Arxiv version.
Funding
I can provide the source code.
## Reattaching the author to the issue ticket: @dl6tud ##
No child items are currently assigned. Use child items to break down this issue into smaller parts.
Link issues together to show that they're related. Learn more.
Activity
- Mark OLESEN assigned to @Sergio
assigned to @Sergio
- Maintainer
Discussed quickly with @andy - is this much different than src/finiteVolume/interpolation/surfaceInterpolation/schemes/harmonic/harmonic.H ?
- Author Maintainer
I denote by phi_1 the potential in cell 1, phi_2 in cell 2, and phi_f the potential at the face between both cells. Similarly, delta_1 means the distance from cell centre 1 to the face, delta_2 the distance of cell 2 to the face, and delta the distance between both cell centres.
Harmonic interpolation is then
φf=φ1δ1/δ+φ2δ2/δ1.However, the missing interpolation scheme is
φf=f⋅φ1+(1−f)⋅φ2with
f = \frac{\delta_2\sigma_1}{\delta_1\sigma_2 + \delta_2\sigma_1}.So yes, both schemes are different. The new schemes accounts for the change of electric conductivity, while harmonic interpolation does not.
By Norbert Weber on 2019-07-31T11:36:58 (imported from GitLab project)
- Maintainer
## Reattaching the author to the issue ticket: @dl6tud ##
- Please register or sign in to reply
- Developer
Hi Norbert,
Is not the harmonic mean the exact solution for the heat flux?. Which is the name of your mean?.
Thanks
Sergio
- Author Maintainer
Hello Sergio,
the harmonic interpolation is perfect for computing the thermal CONDUCTIVITY at the wall:
laplacianSchemes { laplacian(lambda,T) Gauss harmonic uncorrected; }
However, a new scheme is necessary for computing TEMPERATURE at the face:
gradSchemes { grad(T) Gauss volWeighted "lambda"; }
This scheme reads the thermal conductivity from the object registry ("lambda"). Then, it weights the temperature of neighboring cells using the thermal conductivity to get the temperature on the face.
By Norbert Weber on 2019-08-01T15:38:54 (imported from GitLab project)
- Developer
Hi Norbert,
OK. I see your point. We need a correct name for this weighted interpolation. We have a "weighted" interpolation already. Do you know the name of this interpolation? If you provide the code, we can included. Thanks Sergio
- Author Maintainer
I discussed with the colleagues. It is not easy to find a good name for the scheme. The words "weighted" and "flux" (referring to e.g. heatflux) fit well. "Conductivity" does not fit well, because the scheme might be used for diffusivities, as well. Might "weightedFlux" be an option?
By Norbert Weber on 2019-08-07T14:52:14 (imported from GitLab project)
- Author Maintainer
A suggestion for the name:
the scheme you are talking about is required to correctly deal with problems where the field which multiplies the gradient is discontinuous across the face. Otherwise, if this field changes continuously in the domain, the standard linear interpolation would be more appropriate.
Therefore, I’d look for a name which highlights this.
There are also other applications where this scheme can be potentially employed, for example the discretisation of the product 1/rho*grad(p) in the momentum equation when the cell face overlaps with the interface between two fluid phases (see [1]).
[1] Vukčević, V., Jasak, H., Gatin, I., 2017. Implementation of the ghost fluid method for free surface flows in polyhedral finite volume framework. Comput. Fluids 153, 1–19.
By Luca Cornolti on 2019-08-08T12:22:38 (imported from GitLab project)
- Developer
It seems to me that this interpolation scheme is useful for getting fluxes in general such as (flux = K grad(T)), then Tf is needed weighted with K being included. Although the expression fvc::grad(T), k is not specified. It looks like weightedFlux could be an option.
- Author Maintainer
I agree with both of you. For me, "weightedFlux" would be fine. If you agree, I will push the feature at the beginning of September, as I am on holidays for the next two weeks.
By Norbert Weber on 2019-08-09T11:51:38 (imported from GitLab project)
- Maintainer
- Developer
Dear Norbet, On which branch are you trying to push? Can you send me the code to review and I could add it to our development branch Thanks
- Author Maintainer
Dear Sergio,
starting from the development branch, I created a new feature branch and tried to push that (did not work). Please find the code attached to this post. I read the code styling guide, and hope to have respected everything. If I did something wrong, I would be happy for feedback. Thank you!
weightedFluxInterpolationScheme.zip
By Norbert Weber on 2019-09-04T10:25:32 (imported from GitLab project)
- Developer
Dear Norbert, Thanks for the code. It looks good. Still need to check how we want to implement this variation of grad. I think this should look like: fvc::grad(sigma, T) similar as for the fvc::laplacian(kappa, T).
q = sigma*fvc::grad(T), this implies the flux-weighted interpolation of sigma for to get an accurate q, and the user needs to know that the scheme used to interpolate grad(T) uses sigma.
But, ideally we'd want i.e q = fvc::grad(sigma, T), where sigma is already in the expression. where sigma could be a face, cell centre or scalar values or tmp's of those. This requires more coding as the options are larger.
Thanks
Sergio
- Developer
Hi Norbert, I am pushing the code for the new version, I hope. Do you have a simple test case where it can be shown the benefits of this scheme for fluxes? Thanks Sergio
- Author Maintainer
Thank you! Sorry for the silence - I was on holidays. I do have indeed a test case. Should it be for you only, or to push it with the new openfoam version?
By Norbert Weber on 2019-09-28T15:51:27 (imported from GitLab project)
- Developer
I build a new branch with some other numeric. Let me know how you want your contribution in the header file. If you please can send me a small tutorial case which I can add to our tutorial cases. Before making the release we need to follow some code review internally. We'll try to make for this release. Thanks
- Author Maintainer
Dear Sergio, for the best (=simplest) tutorial case I would need laplacianFoam. Unfortunately, laplacianFoam solves
fvm::laplacian(DT, T)
with DT being a dimensionedScalar - but I would need DT as a volScalarField. Might we change laplacianFoam such that it reads DT with READ_IF_PRESENT as volScalarField?It would be good for me if my name, and a reference to the article appears in the header file. Thank you!
By Norbert Weber on 2019-10-03T08:05:00 (imported from GitLab project)
- Developer
Yes, I agree. When you have it ready send it over and I consult if we will make DT a volScalarField instead of dimensionedScalar. Thanks
- Developer
I'd need an affilation for the header:
Copyright (C) YEAR AUTHOR,AFFILIATION
The reference goes bellow.
Thanks
- Author Maintainer
Copyright (C) 2019 Norbert Weber, Helmholtz-Zentrum Dresden-Rossendorf
Please find the testcase here: testcase.zip. It computes the heat flow in a 1D bar, and illustrates the difference of both
laplacian(DT,T) Gauss linear corrected
andlaplacian(DT,T) Gauss harmonic corrected
, andgrad(T) Gauss linear
andgrad(T) Gauss weightedFlux "DT"
. Gnuplot is used to show the results.I modified in laplacianFoam in createFields.H:
dimensionedScalar DTscalar(transportProperties.lookupOrDefault<dimensionedScalar>("DT",0.0)); volScalarField DT ( IOobject ( "DT", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE ), mesh, DTscalar ); volVectorField q ( IOobject ( "q", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE ), mesh, dimensionedVector("q", dimensionSet(0,1,-1,1,0,0,0),vector::zero) );
and in laplacianFoam.C:
q = DT*fvc::grad(T);
Thank you!!!
By Norbert Weber on 2019-10-05T09:28:54 (imported from GitLab project)
- Developer
- Developer
Hi, My mistake, now I can see the flux doing the correct interpolation Thanks
- Developer
- Author Maintainer
Perfect! Thanks a lot!
By Norbert Weber on 2019-10-07T19:26:48 (imported from GitLab project)
- Kutalmış Berçin changed the description
changed the description
- Mark OLESEN mentioned in commit 01b4519f
mentioned in commit 01b4519f
- Mark OLESEN mentioned in commit d04de0d4
mentioned in commit d04de0d4
- Mark OLESEN mentioned in commit 3b6ba059
mentioned in commit 3b6ba059
- Owner
Resolved in v1912 release
- Andrew Heather closed
closed