Just thought I would report back that I've now implemented this way of
copying between systems (even differently sized systems)... and it worked
out to be a nice little function I thought others would appreciate. I use
it twice... once to copy down to a smaller system and once to copy the
solution of the smaller system back out to the larger system.
Enjoy...
void copyVarValues(MeshBase & mesh,
const unsigned int from_system, const unsigned int
from_var, const NumericVector<Number> & from_vector,
const unsigned int to_system, const unsigned int
to_var, NumericVector<Number> & to_vector)
{
MeshBase::node_iterator it = mesh.local_nodes_begin();
MeshBase::node_iterator it_end = mesh.local_nodes_end();
for(;it!=it_end;++it)
{
Node * node = *it;
//The zeroes are for the component.
//If we ever want to use non-lagrange elements we'll have to change
that.
unsigned int from_dof = node->dof_number(from_system,from_var,0);
unsigned int to_dof = node->dof_number(to_system,to_var,0);
to_vector.set(to_dof,from_vector(from_dof));
}
}
Derek
On Wed, Jan 21, 2009 at 3:22 PM, Derek Gaston <[email protected]> wrote:
> On Jan 21, 2009, at 2:56 PM, Roy Stogner wrote:
>
> If all these systems are on the same mesh, no communication should be
>> necessary. You'll still have to use PETSc APIs, but the existing
>> NumericVector::insert() interface might suffice and stay
>> Trilinos-compatible.
>>
>
> Right. They are on the same mesh (all of the systems are part of the same
> EquationSystems object).... and there won't be any parallel communication.
>
>
>>> Not if you do the insert()s element by element. (or
>> element-by-element for elem dofs and node by node for node dofs)
>>
>
> For now... I only have nodal degrees of freedom (no real plans to use
> non-lagrange elements at the current time). So you're saying I can do a
> node loop and call node.dof_number(system, var, 0) (0 since I'm using
> Lagrange) to get the dof_number for that variable in each system.... then
> just do the copy. That's really not bad... and pretty straight forward. I
> suppose I will need to check to make sure that that dof_number is owned by
> the local processor... but that's not a big deal. I guess I can just check
> the processor_id() of the node to see if it corresponds to the local
> processor id before I lookup the dof_numbers...
>
> I'll probably switch to this... what I'm doing does work... but is Petsc
> specific.... and has some fragility built in.
>
>
>>> If I was doing this, I'd first extend System::project_vector() to take
>> a variable specifying projection of only one variable at a time, then
>> if I was worried about efficiency I'd add another overload replacing
>> fptr and gptr with identification of another system and variable so as
>> to avoid the MeshFunction construction and octree lookups. But that
>> would only count as an easier way if it already existed, and "better"
>> depends on whether you like or dislike "more flexible and less
>> efficient".
>>
>
> Hmmm.... I don't know about this. The straight forward method above sounds
> like the way to go.
>
> Thanks!
> Derek
>
------------------------------------------------------------------------------
Create and Deploy Rich Internet Apps outside the browser with Adobe(R)AIR(TM)
software. With Adobe AIR, Ajax developers can use existing skills and code to
build responsive, highly engaging applications that combine the power of local
resources and data with the reach of the web. Download the Adobe AIR SDK and
Ajax docs to start building applications today-http://p.sf.net/sfu/adobe-com
_______________________________________________
Libmesh-users mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/libmesh-users