Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Question in structural/example_3 #27

Open
anupzope opened this issue Aug 15, 2019 · 13 comments
Open

Question in structural/example_3 #27

anupzope opened this issue Aug 15, 2019 · 13 comments

Comments

@anupzope
Copy link

What is the purpose of the loop here?

@manavbhatia
Copy link
Member

This code identifies the node closest to the node pt, which is the location of the center of the panel.
The for loop iterates over all nodes in the mesh and ends with the node closest to pt.

Later this node is used to identify the degree-of-freedom of the the z-displacement at this node, which is then used to collect the value of this dof on rank 0, so that the displacement can be printed at each load-step here.

@anupzope
Copy link
Author

If the mesh object is distributed, after the loop, would this code require a global reduction to get the node with smallest norm w.r.t. the plate's centre?

@manavbhatia manavbhatia reopened this Aug 16, 2019
@manavbhatia
Copy link
Member

If you are using a ReplicatedMesh object then each processor stores the entire mesh and I would expect each rank to compute the same result.

if you are using a DistributedMesh then the result will be different on each processor, and a reduction will be needed to identify the correct node. In this example only rank 0 needs to know about the closest node for output.

@anupzope
Copy link
Author

Do the nonlinear_system and system variables reference the same object?

@manavbhatia
Copy link
Member

manavbhatia commented Aug 16, 2019 via email

@anupzope
Copy link
Author

What is the time step for each iteration? How to tell MAST to use variable time step?

@manavbhatia
Copy link
Member

This example uses a static solve per load step.
For transient solvers you can see the setup in the fluid example.

When setting up the transient solve for structural system you will need to use the second order Newmark solver and the structural transient element assembly objects.

@anupzope
Copy link
Author

How to specify initial velocity of the plate nodes?

@manavbhatia
Copy link
Member

The transient solver stores three vectors for current and previous time-steps: displacement, velocity and acceleration. You can get the velocity vector from the solver API and set the values for the respective degrees-of-freedom. Note that the dof location (that is the location in the vector) is the same for displacement and velocity.

@anupzope
Copy link
Author

Ok. So, I found MAST::TransientSolverBase::solution(), MAST::TransientSolverBase::velocity(), and MAST::TransientSolverBase::acceleration() methods. They return libMesh::NumericVector<Real> reference. However, I am not clear at which index positions are the x, y, z dofs for, let's say i th node, stored. Also, How do these index positions change when using a DistributedMesh?

Since there are three dofs associated with each node, I am assuming that libMesh::NumericVector<Real>::last_local_index()-libMesh::NumericVector<Real>::first_local_index() is three times the number of local nodes on an MPI process. Please correct me if I am wrong.

@manavbhatia
Copy link
Member

Ok. So, I found MAST::TransientSolverBase::solution(), MAST::TransientSolverBase::velocity(), and MAST::TransientSolverBase::acceleration() methods. They return libMesh::NumericVector<Real> reference.

Correct.

However, I am not clear at which index positions are the x, y, z dofs for, let's say i th node, stored. Also, How do these index positions change when using a DistributedMesh?

This comes from the node pointer. libMesh::Node inherits from libMesh::DofObject, which stores the dof information. You can call :

  • libMesh::Node::dof_number(0, 0, 0) to get the index for u on the node,
  • libMesh::Node::dof_number(0, 1, 0) to get the index for v on the node,
  • libMesh::Node::dof_number(0, 2, 0) to get the index for w on the node.

These identify the location of the corresponding value for the node in the NumericVector.
In a parallel run each rank will set these values based on the first and last local index of the vector.

Since there are three dofs associated with each node,

Actually, there are 6 dofs on each node: u,v,w, tx, ty, tz.

I am assuming that libMesh::NumericVector<Real>::last_local_index()-libMesh::NumericVector<Real>::first_local_index() is three times the number of local nodes on an MPI process. Please correct me if I am wrong.

This will be 6 times number of local nodes.

@anupzope
Copy link
Author

Thanks for the correction.

So, is following code correct to set zero initial displacement and acceleration, but non-zero initial velocity?

typedef libMesh::MeshBase::const_node_iterator nodeiter;

libMesh::NumericVector<Real> & disp = solver.solution();
libMesh::NumericVector<Real> & vel = solver.velocity();
libMesh::NumericVector<Real> & acc = solver.acceleration();

disp.zero();
acc.zero();

nodeiter niter = mesh.local_nodes_begin();
nodeiter niterend = mesh.local_nodes_end();
Real initialVelocity[3] = {0.0, 0.01, 0.0};
for(; niter != niterend; ++niter) {
  libMesh::dof_id_type i = (*niter)->dof_number(0, 0, 0);
  libMesh::dof_id_type j = (*niter)->dof_number(0, 1, 0);
  libMesh::dof_id_type k = (*niter)->dof_number(0, 2, 0);
  
  vel.set(i, initialVelocity[0]);
  vel.set(j, initialVelocity[1]);
  vel.set(k, initialVelocity[2]);
}

@manavbhatia
Copy link
Member

I think this should work. Remember to call NumericVector::close() on all vectors after you are done setting the values.

After setting the initial conditions, you will need to call TransientSolverBase::solve_highest_derivative_and_advance_time_step() see here

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
None yet
Projects
None yet
Development

No branches or pull requests

2 participants