[Nek5000-users] weak derivative, dssum and binvm1 in cvode_driver
nek5000-users at lists.mcs.anl.gov
nek5000-users at lists.mcs.anl.gov
Thu Jun 5 10:11:12 CDT 2014
I guess you want to solve for species and temperature using detailed
chemistry and transport. Right?
Von: <nek5000-users at lists.mcs.anl.gov>
Antworten an: <nek5000-users at lists.mcs.anl.gov>
Datum: Thursday, June 5, 2014 at 5:06 PM
An: "nek5000-users at lists.mcs.anl.gov" <nek5000-users at lists.mcs.anl.gov>
Betreff: Re: [Nek5000-users] weak derivative, dssum and binvm1 in
cvode_driver
Hi Paul and Stefan,
Thank you for the explanations, things are going gradually more clear in my
mind.
So for a practical application, still in the context of cvode, I would like
to add a new term to the RHS of the temperature equation. This new term is
roughly
(vdiff * grad (Y)) grad (T)
where Y is a passive scalar, vdiff the diffusion coefficient of this passive
scalar and T the temperature.
Basically in the routine makeq.f the implementation should looks like to
call wgradm1(Y_x, Y_y, Y_z, t(1,1,1,1,ifield+1), nelv)
call dssum(grad_Y, Y_x, Y_y, Y_z)
call col2(grad_Y,bintm1,ntot)
call col2(grad_Y, vdiff(1,1,1,1,ifield+2) ,ntot )
call wgradm1(T_x, T_y, T_z, t(1,1,1,1,ifield), nelv)
call dssum(grad_T, T_x, T_y, T_z)
call col2(grad_T,grad_Y,ntot)
call add2(bq(1,1,1,1,ifield), grad_T,ntot)
and then convab and wlaplacian routines will add the other terms, and it
will finish with
call dssum(bq,nx1,ny1,nz1)
call col2(bq,bintm1,ntot)
call col2(bq,bm1,ntot)
Am I right??
Thank you for your help !
Best regards,
Emmanuel
On 5 Jun 2014, at 19:25, <nek5000-users at lists.mcs.anl.gov>
<nek5000-users at lists.mcs.anl.gov> wrote:
>
> Hi Emmanuel,
>
> Note that bintm1 and bm1 are not precise inverses of one another.
> bm1 is the unassembled mass matrix - i.e., it is the collection of local
> (to each element) mass matrices, each of which is diagonal.
>
> Conceptually, we denote the local mass matrix as B_L, for local.
> Then, the assembled mass matrix is:
>
> B := Q^T B_L Q
>
> where Q is the Boolean matrix that copies global degrees of freedom
> to their shared counterparts (i.e., copies a values to shared faces,
> edges, vertices, etc.). For the SEM, B_L is diagonal, as is B, so
> Binv := B^-1 is easily computed. Once the entries of this diagonal
> matrix are known, we can copy them to their element-based locations,
> so that we can generate B^-1 * r very easily.
>
> Often, products with these matrices, in conjunction with dssum, are used
> to generate continuous fields and that appears to be what is happening in
> this case. An example would be if you want to find a the continuous
> counterpart to a discontinuous function g. The usual projection approach
> is: Find gc \in H1 such that \int v gc dV = \int v g dV for all v in H1,
> which
> ends up looking like:
>
> Find gc \in H1 such that:
>
>
> (v,gc) = (v,g) for all g \in H1
>
> or
>
> v B gc = v Q^T B_L g_L
>
> gc = B^-1 Q^T B_L g_L
>
> gc_L = Q B^-1 Q^T B_L g_L
>
> In Nek, the last statement is expressed as:
>
> call col2(g,bm1,n)
> call dssum(g,nx1,ny1,nz1)
> call col2(g,bintm1,n)
>
> In your cvode case, these statements are permuted because of the
> initial state of the data and the desired output (which is to be B*q --
> which is why the variable is labeled bq).
>
> Paul
>
>
>
>
>
> From: nek5000-users-bounces at lists.mcs.anl.gov
> [nek5000-users-bounces at lists.mcs.anl.gov] on behalf of
> nek5000-users at lists.mcs.anl.gov [nek5000-users at lists.mcs.anl.gov]
> Sent: Thursday, June 05, 2014 2:38 AM
> To: nek5000-users at lists.mcs.anl.gov
> Subject: Re: [Nek5000-users] weak derivative, dssum and binvm1 in cvode_driver
>
> Hi Emmanuel
>
> From what I recall:
> CVODE solves the governing non-linear eqns in the strong form. That¹s why we
> have to project every term of the RHS onto H2 (dssum and collocate with
> bintm1). The sum is projected again onto H2 because it gets devided by vtrans.
> Hope this helps.
>
> Cheers,
> Stefan
>
> Von: <nek5000-users at lists.mcs.anl.gov>
> Antworten an: <nek5000-users at lists.mcs.anl.gov>
> Datum: Thursday, June 5, 2014 at 9:20 AM
> An: "nek5000-users at lists.mcs.anl.gov" <nek5000-users at lists.mcs.anl.gov>
> Betreff: Re: [Nek5000-users] weak derivative, dssum and binvm1 in cvode_driver
>
> I¹m adding another quick question:
>
> why in the routine makeq, under the cvode statement, there are these following
> 2 lines:
> call col2(bq,bintm1,ntot)
> call col2(bq,bm1,ntot)
>
> I¹m not sure to understand why we are multiplying bq by the inverse of the
> mass matrix (bintm1) and just after multiplying it back by the mass matrix
> (bm1). Moreover, after leaving the makeq routine, fctfun will multiplying
> again ydot=bq/vtrans by the inverse of the mass matrix.
>
> Why all of these operations? Thank you for your help !
>
>
>
> Best regards,
> Emmanuel
>
>
>
> On 5 Jun 2014, at 12:14, <nek5000-users at lists.mcs.anl.gov>
> <nek5000-users at lists.mcs.anl.gov> wrote:
>
>> Dear Nek team,
>>
>> I have a quick question: in cvode_driver.f, in the routine fcvfun, do all the
>> derivatives must be in weak formulation (using wgradm1) when constructing the
>> right hand side term ydot?
>> Because at the end of the routine fcvfun, dssum and multiplication by binvm1
>> are applied to the right hand side term ydot.
>>
>> It¹s also not clear to me if I need to apply dssum and binvm1 for each
>> derivative term appearing in my equation (in other terms, each time after
>> calling wgradm1), or if I need to apply dssum and binvm1 only to the whole
>> RHS term ydot ?
>>
>> I¹m afraid of inconsistency if apply twice dssum and binvm1: one during the
>> construction of operators and one on the final ydot term.
>>
>> Thank you for you help !
>>
>> Best regards,
>> Emmanuel
>>
>>
>>
>> _______________________________________________
>> Nek5000-users mailing list
>> Nek5000-users at lists.mcs.anl.gov
>> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
>
> _______________________________________________ Nek5000-users mailing
> listNek5000-users at lists.mcs.anl.gov
> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
> _______________________________________________
> Nek5000-users mailing list
> Nek5000-users at lists.mcs.anl.gov
> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
_______________________________________________ Nek5000-users mailing list
Nek5000-users at lists.mcs.anl.gov
https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/nek5000-users/attachments/20140605/656929e3/attachment-0001.html>
More information about the Nek5000-users
mailing list