<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=Windows-1252">
</head>
<body style="word-wrap: break-word; -webkit-nbsp-mode: space; -webkit-line-break: after-white-space;">
Hi Paul
<div><br>
</div>
<div>I have some difficulties to compute a simple gradient of a variable. I think I’m confusing something.</div>
<div><br>
</div>
<div>Basically I would like to modify the temperature and species equations to add a new term in order to take into account a velocity correction related to mass conservation and species diffusion coefficients.</div>
<div><br>
</div>
<div>For temperature the equation is:</div>
<div><br>
</div>
<div>d T /dt = term1 + term2 + term3</div>
<div>with</div>
<div>term1 = - U grad (T)</div>
<div>term2 = 1/rho Cp * grad ( lambda * grad (T))</div>
<div>term3 = 1/rho Cp * (rho * sum_k [Cp,k * (Y_k * D_k * grad(Y_k) - Y_k * sum_j [D_j * grad(Y_j )] ) ] ) * grad(T)</div>
<div><br>
</div>
<div>and for species equation :</div>
<div><br>
</div>
<div>
<div>d Y /dt = term1 + term2 + term3</div>
<div>with</div>
<div>term1 = - U grad (Y)</div>
<div>term2 = 1/rho Cp * grad ( lambda * grad (Y))</div>
</div>
<div>term3 = - 1/rho * grad (Y * sum_k [D_k * grad (Y_k) ] )</div>
<div><br>
</div>
<div>where D_k or D_j are the species diffusion coefficients, T the temperature and Y the species mass fraction.</div>
<div><br>
</div>
<div>I am using the low-Mach number formulation and solving the equations with CVODE.</div>
<div>Basically term1 is computed by the routine convab while term2 is computed by routine wlapacian. This is the current official nek implementation.</div>
<div><br>
</div>
<div>Now I would like to add the term3, which involves some gradients of temperature or species variables. I think the simplest way is to modify the routine convab or maybe makeq.f</div>
<div><br>
</div>
<div>However I tried a simple call gradm1 (Tx,Ty,Tz,T) but NEK crashes with a segmentation fault error. I tried to debug the routine gradm1, and it appears that it crashes during the loop </div>
<div><br>
</div>
<div>
<div>do i=1,lxyz</div>
<div> ux(i,e) =jacmi(i,e)*(ur(i)*rxm1(i,1,1,e)</div>
<div> $ + us(i)*sxm1(i,1,1,e) )</div>
<div> uy(i,e) =jacmi(i,e)*(ur(i)*rym1(i,1,1,e)</div>
<div> $ + us(i)*sym1(i,1,1,e) )</div>
<div> enddo</div>
</div>
<div><br>
</div>
<div>after reaching a particular element e (in my case it crashes for element 9).</div>
<div><br>
</div>
<div><br>
</div>
<div>So my questions are just:</div>
<div>1) According to the implementation of term3 in temperature and species equations, what is the best way to proceed? I have think about modifying convab or makeq but is it a good idea?</div>
<div><br>
</div>
<div>2) Is it a good idea to use gradm1? Why does it crash for a particular element? Should I pay attention to something in particular when using this routine?</div>
<div><br>
</div>
<div><br>
</div>
<div>Thank you very much for your help !</div>
<div><br>
</div>
<div><br>
</div>
<div><br>
<div>
<div style="color: rgb(0, 0, 0); letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-stroke-width: 0px; word-wrap: break-word; -webkit-nbsp-mode: space; -webkit-line-break: after-white-space;">
<div>Best regards,</div>
<div>Emmanuel</div>
</div>
<br class="Apple-interchange-newline">
<br class="Apple-interchange-newline">
</div>
<br>
<div>
<div>On 30 May 2014, at 10:25, <a href="mailto:nek5000-users@lists.mcs.anl.gov">nek5000-users@lists.mcs.anl.gov</a> wrote:</div>
<br class="Apple-interchange-newline">
<blockquote type="cite">
<div style="font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; line-height: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-stroke-width: 0px; direction: ltr; font-family: Tahoma; font-size: 10pt;">
<br class="Apple-interchange-newline">
Hi Emmanuel,
<div><br>
</div>
<div>Strictly speaking, our functions are only once differentiable -- they are C0, not C1.</div>
<div><br>
</div>
<div>The test functions are also in H1 (i.e., are C0) and therefor the standard Galerkin</div>
<div>approach is to transfer one derivative to the test functions via integration by parts</div>
<div>and thus avoid twice differentiating a function.</div>
<div><br>
</div>
<div>If you just want the gradient of something, you would use gradm1 (gradient for</div>
<div>a quantity on mesh 1). </div>
<div><br>
</div>
<div>Paul</div>
<div><br>
</div>
<div><br>
<div style="font-family: 'Times New Roman'; font-size: 16px;">
<hr tabindex="-1">
<div id="divRpF290929" style="direction: ltr;"><font face="Tahoma" size="2"><b>From:</b><span class="Apple-converted-space"> </span><a href="mailto:nek5000-users-bounces@lists.mcs.anl.gov">nek5000-users-bounces@lists.mcs.anl.gov</a><span class="Apple-converted-space"> </span>[<a href="mailto:nek5000-users-bounces@lists.mcs.anl.gov">nek5000-users-bounces@lists.mcs.anl.gov</a>]
on behalf of<span class="Apple-converted-space"> </span><a href="mailto:nek5000-users@lists.mcs.anl.gov">nek5000-users@lists.mcs.anl.gov</a><span class="Apple-converted-space"> </span>[<a href="mailto:nek5000-users@lists.mcs.anl.gov">nek5000-users@lists.mcs.anl.gov</a>]<br>
<b>Sent:</b><span class="Apple-converted-space"> </span>Thursday, May 29, 2014 7:47 PM<br>
<b>To:</b><span class="Apple-converted-space"> </span><a href="mailto:nek5000-users@lists.mcs.anl.gov">nek5000-users@lists.mcs.anl.gov</a><br>
<b>Subject:</b><span class="Apple-converted-space"> </span>Re: [Nek5000-users] Questions in AXHELM and WLAPLACIAN<br>
</font><br>
</div>
<div></div>
<div>Hi Paul,
<div><br>
</div>
<div>Thank you for your quick answer, this is really helpful.</div>
<div><br>
</div>
<div>I still have a question concerning the computation of a gradient.</div>
<div><br>
</div>
<div>Why did you construct the diffusion operator with the helmholtz routine? Why not using directly the routines gradm1 ?</div>
<div><br>
</div>
<div><br>
</div>
<div>By the way I am not sure to understand the differences between gradm1 (navier5.f), wgradm1 (navier1.f) , gradm11 (navier2.f) and gradm11ts (navier2.f).</div>
<div>I guess that the w letter refers to the weak formulation, but I don’t understand in which cases we need to choose a formulation instead of others.</div>
<div><br>
</div>
<div>In summary, if I want to modify the temperature/passive scalar equations to add a term that contains a gradient of the solved temperature/PS variable, for example, which “gradm” routine do I need to use?</div>
<div><br>
</div>
<div>Thank you for your help !</div>
<div><br>
</div>
<div><br>
<div>
<div style="letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; word-wrap: break-word;">
<div>Best regards,</div>
<div>Emmanuel</div>
</div>
<br class="Apple-interchange-newline">
<br class="Apple-interchange-newline">
</div>
<br>
<div>
<div>On 29 May 2014, at 20:20,<span class="Apple-converted-space"> </span><a href="mailto:nek5000-users@lists.mcs.anl.gov" target="_blank">nek5000-users@lists.mcs.anl.gov</a><span class="Apple-converted-space"> </span>wrote:</div>
<br class="Apple-interchange-newline">
<blockquote type="cite"><br>
<br>
Hi Emmanuel,<br>
<br>
axhelm effects the weak Laplacian plus the mass matrix.<br>
<br>
Specifically, if you wanted to solve:<br>
<br>
(1) -nabla . [ h1 * nabla ] u + h2 * u = f<br>
<br>
with h1, h2, and f all functions of X=(x,y,z),<br>
the discrete form would be:<br>
<br>
<br>
(2) A(h1) u + B(h2) u = B f<br>
<br>
or<br>
<br>
(3) H u = B f<br>
<br>
where A is the (h1-dependent stiffness matrix and<br>
B(h2) is the h2-dependent mass matrix.<br>
<br>
For the SEM, B is diagonal, so B(h2) = h2*B, where<br>
we view h2 as a diagonal matrix.<br>
<br>
In Nek5000, a matrix-vector product w = Hu is effected via<br>
<br>
call axhelm (w,u,h1,h2,imesh,isd)<br>
<br>
with w,u,h1,h2 being arrarys. Here, imesh indicates<br>
if you are on the fluid mesh or thermal mesh. The fluid<br>
mesh is a subset of the thermal mesh in the case of conjugate heat transfer between a fluid and a solid.<br>
isd indicates which spatial direction is considered<br>
in axisymmetric applications where u is a component<br>
of a vector field.<br>
<br>
Notice the sign of the Laplacian in Eq. (1) -- that may<br>
help resolve your second question? (There is a sign<br>
change that ariss from the standard integration-by-parts<br>
stepin the weak form.)<br>
<br>
Paul<br>
<br>
<br>
On Thu, 29 May 2014,<span class="Apple-converted-space"> </span><a href="mailto:nek5000-users@lists.mcs.anl.gov" target="_blank">nek5000-users@lists.mcs.anl.gov</a><span class="Apple-converted-space"> </span>wrote:<br>
<br>
<blockquote type="cite">Dear nek5000 team,<br>
<br>
Could you help me with the 2 questions below? Thank you very much !<br>
<br>
1) Concerning the routine “axhelm” in hmholtz.f. If I am right, one<br>
possibility of this routine can be the discretisation of the diffusion<br>
operator. However it is not clear to me if axhelm performs:<br>
<br>
a) vdiff x Laplacian(phi)<br>
<br>
or<br>
<br>
b) grad ( vdiff x grad (phi) )<br>
<br>
phi being the variable solved for (T or PS for example).<br>
<br>
So, is it actually a) or b) ?<br>
<br>
<br>
2) In the routine wlaplacian in navier1.f, I am not sure to understand the<br>
consequence of the line:<br>
<br>
call sub2 (out,wrk,ntot)<br>
<br>
after the lines:<br>
<br>
call bcneusc(out,1) call axhelm(wrk,a,diff,h2,imesh,1)<br>
<br>
If I am right, this will change the sign in front of the diffusion operator,<br>
right? Because in makeq.f, I am not sure that the line call<br>
add2(bq(1,1,1,1,ifield-1),w1,ntot)<br>
<br>
will actually assemble the RHS of the transport equation : - (vtrans) x U<br>
grad (phi) + grad ( vdiff x grad (phi) )<br>
<br>
<br>
Thank you for your help !<br>
<br>
Best regards, Emmanuel<br>
<br>
</blockquote>
_______________________________________________<br>
Nek5000-users mailing list<br>
<a href="mailto:Nek5000-users@lists.mcs.anl.gov" target="_blank">Nek5000-users@lists.mcs.anl.gov</a><br>
<a href="https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users">https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users</a><br>
</blockquote>
</div>
<br>
</div>
</div>
</div>
</div>
</div>
<span style="font-family: Helvetica; font-size: 12px; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; line-height: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-stroke-width: 0px; float: none; display: inline !important;">_______________________________________________</span><br style="font-family: Helvetica; font-size: 12px; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; line-height: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-stroke-width: 0px;">
<span style="font-family: Helvetica; font-size: 12px; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; line-height: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-stroke-width: 0px; float: none; display: inline !important;">Nek5000-users
mailing list</span><br style="font-family: Helvetica; font-size: 12px; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; line-height: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-stroke-width: 0px;">
<a href="mailto:Nek5000-users@lists.mcs.anl.gov" style="font-family: Helvetica; font-size: 12px; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; line-height: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-stroke-width: 0px;">Nek5000-users@lists.mcs.anl.gov</a><br style="font-family: Helvetica; font-size: 12px; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; line-height: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-stroke-width: 0px;">
<a href="https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users" style="font-family: Helvetica; font-size: 12px; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; line-height: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-stroke-width: 0px;">https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users</a></blockquote>
</div>
<br>
</div>
</body>
</html>