Dear Birger,

thanks for this comment, you are right.

You can directly use *scv.volume()* to get the volume.

Regards,
Martin


On 08/22/2018 10:38 AM, Birger Hagemann wrote:

Dear Martin,

The unit for the well index is m³. Consequently the unit for your expression for values[contiWEEqIdx] is kg/s. However, the source term has the unit kg/m³/s. I think you have to divide by the cell volume:

*values[contiWEqIdx] = densityW/viscosityW * WI *(pbh - pw)/fvGeometry.subContVol[scvIdx].volume;*

**

Am I wrong?

Regards Birger

*Von:*Dumux [mailto:[email protected]] *Im Auftrag von *Martin Schneider
*Gesendet:* Mittwoch, 22. August 2018 10:03
*An:* DuMuX User Forum <[email protected]>; Ranjeet Singh <[email protected]>
*Betreff:* Re: [DuMuX] How to implement Peaceman well model in DuMuX?

Dear Ranjeet,

you can implement the Peaceman well model as a solution-dependent source,
using the following function in your problem file

*NumEqVector source(const Element &element,
                                       const FVElementGeometry& fvGeometry,                                        const ElementVolumeVariables& elemVolVars,
                                       const SubControlVolume &scv) const*

(see for example porousmediumflow/2pnc/implicit/fuelcellproblem.hh)

With the *elemVolVars *you have access to the pressure, density, etc.:
*const auto& volVars = elemVolVars[scv];
Scalar pw = volVars.pressure(wPhaseIdx);
Scalar densityW = volVars.density(wPhaseIdx);
Scalar viscosityW = volVars.viscosity(wPhaseIdx);
*
Such that the source term looks as follows (if you inject water and neglect gravity)
*NumEqVector values(0.0);*
*values[contiWEqIdx] = densityW/viscosityW * WI *(pbh - pw);

*where *WI *is the well index (see the Peaceman well model); *pbh* the bottom hole pressure, etc.
If gravity is included you have to use the phase potentials.

Kind regards,
Martin*

*
On 08/09/2018 04:03 AM, Ranjeet Singh wrote:

    Hi,

      I want use Peaceman well model in DuMuX?.  Could you please some
    information (or any example if already implemented) so that I can
    implement the well model?

    Thank you!

    Regards,

    Ranjeet




    _______________________________________________

    Dumux mailing list

    [email protected]
    <mailto:[email protected]>

    https://listserv.uni-stuttgart.de/mailman/listinfo/dumux

--
M.Sc. Martin Schneider
University of Stuttgart
Institute for Modelling Hydraulic and Environmental Systems
Department of Hydromechanics and Modelling of Hydrosystems
Pfaffenwaldring 61
D-70569 Stuttgart
Tel: (+49) 0711/ 685-69159
Fax: (+49) 0711/ 685-60430
E-Mail:[email protected] <mailto:[email protected]>


_______________________________________________
Dumux mailing list
[email protected]
https://listserv.uni-stuttgart.de/mailman/listinfo/dumux


--
M.Sc. Martin Schneider
University of Stuttgart
Institute for Modelling Hydraulic and Environmental Systems
Department of Hydromechanics and Modelling of Hydrosystems
Pfaffenwaldring 61
D-70569 Stuttgart
Tel: (+49) 0711/ 685-69159
Fax: (+49) 0711/ 685-60430
E-Mail: [email protected]

_______________________________________________
Dumux mailing list
[email protected]
https://listserv.uni-stuttgart.de/mailman/listinfo/dumux

Reply via email to