<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=UTF-8">
</head>
<body text="#000000" bgcolor="#FFFFFF">
<div class="moz-cite-prefix">Hi Seongmo,</div>
<div class="moz-cite-prefix"><br>
</div>
<div class="moz-cite-prefix">I am using following steps to add
annotation to Catalyst,</div>
<div class="moz-cite-prefix"><br>
</div>
<div class="moz-cite-prefix">1 - Create a ProgrammableFilter<br>
2 - Add following code to it (don’t forget to modify epoc)<br>
<br>
import paraview.simple<br>
import paraview.servermanager<br>
import datetime<br>
import locale<br>
<br>
# get output port and copy it to input<br>
out = self.GetOutput()<br>
out.ShallowCopy(self.GetInput())<br>
<br>
# define epoc<br>
locale.setlocale(locale.LC_TIME, "en_US")<br>
epoc = datetime.datetime(1949, 12, 1, 0, 0, 0)<br>
<br>
# get time value<br>
tval = self.GetInput().GetInformation().Get(out.DATA_TIME_STEP())<br>
<br>
# find time by adding hours<br>
#t1 = epoc+datetime.timedelta(seconds=tval)<br>
t1 = epoc+datetime.timedelta(hours=tval)<br>
tstr = t1.strftime(“%d-%b-%Y %H:%M:%S UTC")<br>
<br>
# output<br>
label = vtk.vtkStringArray()<br>
label.SetNumberOfComponents(1)<br>
label.Resize(1)<br>
label.SetName("TimeLabel")<br>
label.SetValue(0, t1)<br>
out.GetFieldData().AddArray(label)<br>
<br>
3 - Add AnnotateGlobalData1 and select “TimeLabel†as array (needs
MergeBlocks filter in case of multiblock datasets)</div>
<div class="moz-cite-prefix"><br>
</div>
<div class="moz-cite-prefix">You could modify the section related
with generation of time string for your needs. I hope it helps.</div>
<div class="moz-cite-prefix"><br>
</div>
<div class="moz-cite-prefix">Regards,</div>
<div class="moz-cite-prefix"><br>
</div>
<div class="moz-cite-prefix">--ufuk<br>
</div>
<div class="moz-cite-prefix"><br>
</div>
<div class="moz-cite-prefix">On 25.01.2018 00:56, Andy Bauer wrote:<br>
</div>
<blockquote type="cite"
cite="mid:CAMaOp+FqVGEEEwQxywUDqD0iR+m0HJKUokedO3vLuzW3BAKNfQ@mail.gmail.com">
<div dir="ltr">
<div>
<div>
<div>Hi Seongmo,<br>
<br>
</div>
Please keep the conversations on the mailing list so that
anyone can follow along or participate. Also, these types of
things often get lost in my inbox when they don't make it
back to the ParaView mailing list.<br>
<br>
</div>
What version of ParaView Catalyst are you using? I think the
annotate time filter should work with Catalyst but I haven't
verified that. I vaguely remember others using that filter
with Catalyst though. Also, I think the colormap bug was
fixed. If you have a way of sharing a sample that demonstrates
either of those bugs I can try taking a look at the issue.<br>
<br>
</div>
<div>Best,<br>
</div>
<div>Andy <br>
</div>
<div>
<div>
<div class="gmail_extra"><br>
<div class="gmail_quote">On Thu, Jan 18, 2018 at 6:40 PM,
Seong Mo Yeon <span dir="ltr"><<a
href="mailto:seongmo.yeon@gmail.com" target="_blank"
moz-do-not-send="true">seongmo.yeon@gmail.com</a>></span>
wrote:<br>
<blockquote class="gmail_quote" style="margin:0 0 0
.8ex;border-left:1px #ccc solid;padding-left:1ex">
<div>
<div name="messageBodySection">Dear Andy Bauer<br>
<br>
I have a quick question.<br>
Is it possible to have annotate time filter
processed in catalyst adaptor? Current my code
cannot that filter.<br>
<br>
BTW, image extracted from catalyst looks different
from render view of paraview at the time of
writing a script. e.g., pressure colormap legend
is missing.<br>
<br>
Regards<br>
Seongmo</div>
<div name="messageReplySection"><br>
On 2018ë…„ 1ì›” 18ì¼ AM 1:17 +0900, Andy Bauer <<a
href="mailto:andy.bauer@kitware.com"
target="_blank" moz-do-not-send="true">andy.bauer@kitware.com</a>>,
wrote:<br>
<blockquote type="cite">
<div dir="ltr">
<div>Hi,<br>
<br>
</div>
My guess is that the TimeStep isn't getting
set properly in the adaptor (though it looks
like it should be in
"dataDescription->SetTimeData(<wbr>runTime.value(),
runTime.deltaTValue());"). My suggestion would
be to add in the following to either the
RequestDataDescription() or DoCoProcessing()
methods in the python script to see what
Catalyst thinks the time step is:<br>
   print("In script2.py, the data time step
is ", datadescription.GetTimeStep())<br>
<br>
</div>
<div class="gmail_extra"><br>
<div class="gmail_quote">On Wed, Jan 17, 2018
at 9:57 AM, SeongMo <span dir="ltr"><<a
href="mailto:seongmo.yeon@gmail.com"
target="_blank" moz-do-not-send="true">seongmo.yeon@gmail.com</a>></span>
wrote:<br>
<blockquote class="gmail_quote"
style="margin:0 0 0 .8ex;border-left:1px
#ccc solid;padding-left:1ex">
<div class="m_-8627969882828805903HOEnZb">
<div class="m_-8627969882828805903h5">Hi,<br>
<br>
I wrote a OpenFOAM adaptor for
Catalyst.<br>
<br>
In the ParaView, the connection is
made good and shows filtered flow
field as written in the python script.<br>
<br>
However, filename_%t and image_%t is
not expanded as time marching but just
write filename_0 and image_0.png.<br>
<br>
As far as I know, %t should be
replaced with current time as given in
dataDescription->SetTimeData.<br>
<br>
Any help would be appreciated.<br>
<br>
<br>
FYI, python script is attached and
snippet of my OpenFOAM Adaptor code
for Catalyst is as follows:<br>
<br>
// icoFoam.C<br>
<br>
#ifdef USE_CATALYST<br>
   Foam::HashTable<string>
options = args.options();<br>
   IStringStream
is(options["scriptList"]);<br>
   wordList scriptList =
readList<word>(is);<br>
   OFAdaptor::Initialize(scriptLi<wbr>st,
mesh);<br>
#endif<br>
   while (runTime.loop())<br>
   {<br>
       runTime.write();<br>
#ifdef USE_CATALYST<br>
       OFAdaptor::CoProcess(mesh,
runTime);<br>
#endif<br>
   }<br>
#ifdef USE_CATALYST<br>
   OFAdaptor::Finalize();<br>
#endif<br>
<br>
<br>
// OFAdaptor.C<br>
<br>
void CoProcess(Foam::fvMesh& mesh,
Foam::Time& runTime)<br>
{<br>
     Â
vtkNew<vtkCPDataDescription>
dataDescription;<br>
     Â
dataDescription->AddInput("inp<wbr>ut");<br>
     Â
dataDescription->SetTimeData(r<wbr>unTime.value(),
runTime.deltaTValue());<br>
      if (runTime.end())<br>
      {<br>
          // assume that we want to
all the pipelines to execute<br>
          // if it is the last time
step<br>
         Â
dataDescription->ForceOutputOn<wbr>();<br>
      }<br>
      if
(Processor->RequestDataDescrip<wbr>tion(dataDescription.GetPointe<wbr>r())
!= 0)<br>
      {<br>
         Â
Foam::polyMesh::readUpdateStat<wbr>e
meshState = mesh.readUpdate();<br>
<br>
          if(meshState !=
Foam::polyMesh::UNCHANGED)<br>
          {<br>
              BuildVTKGrid(mesh);<br>
          }<br>
          UpdateVTKAttributes(mesh);<br>
dataDescription->GetInputDescr<wbr>iptionByName("input")->SetGrid<wbr>(multiBlockDataSet);<br>
         Â
Processor->CoProcess(dataDescr<wbr>iption.GetPointer());<br>
      }<br>
}<br>
<br>
<br>
--<br>
SeongMo Yeon, Ph.D, Senior Engineer<br>
Offshore Hydrodynamics Research<br>
SAMSUNG HEAVY INDUSTRIES CO., LTD.<br>
Central Research Institute<br>
E-mail : <a
href="mailto:seongmo.yeon@gmail.com"
target="_blank"
moz-do-not-send="true">seongmo.yeon@gmail.com</a><br>
Tel :<br>
------------------------------<wbr>--------------------------<br>
Fluctuat nec mergitur<br>
<br>
</div>
</div>
<br>
______________________________<wbr>_________________<br>
Powered by <a
href="http://www.kitware.com"
rel="noreferrer" target="_blank"
moz-do-not-send="true">www.kitware.com</a><br>
<br>
Visit other Kitware open-source projects
at <a
href="http://www.kitware.com/opensource/opensource.html"
rel="noreferrer" target="_blank"
moz-do-not-send="true">http://www.kitware.com/opensou<wbr>rce/opensource.html</a><br>
<br>
Please keep messages on-topic and check
the ParaView Wiki at: <a
href="http://paraview.org/Wiki/ParaView"
rel="noreferrer" target="_blank"
moz-do-not-send="true">http://paraview.org/Wiki/ParaV<wbr>iew</a><br>
<br>
Search the list archives at: <a
href="http://markmail.org/search/?q=ParaView"
rel="noreferrer" target="_blank"
moz-do-not-send="true">http://markmail.org/search/?q=<wbr>ParaView</a><br>
<br>
Follow this link to subscribe/unsubscribe:<br>
<a
href="https://paraview.org/mailman/listinfo/paraview"
rel="noreferrer" target="_blank"
moz-do-not-send="true">https://paraview.org/mailman/l<wbr>istinfo/paraview</a><br>
<br>
</blockquote>
</div>
<br>
</div>
</blockquote>
</div>
</div>
</blockquote>
</div>
<br>
</div>
</div>
</div>
</div>
<br>
<fieldset class="mimeAttachmentHeader"></fieldset>
<pre class="moz-quote-pre" wrap="">_______________________________________________
Powered by <a class="moz-txt-link-abbreviated" href="http://www.kitware.com">www.kitware.com</a>
Visit other Kitware open-source projects at <a class="moz-txt-link-freetext" href="http://www.kitware.com/opensource/opensource.html">http://www.kitware.com/opensource/opensource.html</a>
Please keep messages on-topic and check the ParaView Wiki at: <a class="moz-txt-link-freetext" href="http://paraview.org/Wiki/ParaView">http://paraview.org/Wiki/ParaView</a>
Search the list archives at: <a class="moz-txt-link-freetext" href="http://markmail.org/search/?q=ParaView">http://markmail.org/search/?q=ParaView</a>
Follow this link to subscribe/unsubscribe:
<a class="moz-txt-link-freetext" href="https://paraview.org/mailman/listinfo/paraview">https://paraview.org/mailman/listinfo/paraview</a>
</pre>
</blockquote>
<p><br>
</p>
</body>
</html>